1 #include <petscdmshell.h> /*I "petscdmshell.h" I*/ 2 #include <petscmat.h> 3 #include <petsc/private/dmimpl.h> 4 5 typedef struct { 6 Vec Xglobal; 7 Vec Xlocal; 8 Mat A; 9 VecScatter gtol; 10 VecScatter ltog; 11 VecScatter ltol; 12 void *ctx; 13 } DM_Shell; 14 15 /*@ 16 DMGlobalToLocalBeginDefaultShell - Uses the GlobalToLocal VecScatter context set by the user to begin a global to local scatter 17 Collective 18 19 Input Arguments: 20 + dm - shell DM 21 . g - global vector 22 . mode - InsertMode 23 - l - local vector 24 25 Level: advanced 26 27 Note: This is not normally called directly by user code, generally user code calls DMGlobalToLocalBegin() and DMGlobalToLocalEnd(). If the user provides their own custom routines to DMShellSetLocalToGlobal() then those routines might have reason to call this function. 28 29 .seealso: DMGlobalToLocalEndDefaultShell() 30 @*/ 31 PetscErrorCode DMGlobalToLocalBeginDefaultShell(DM dm,Vec g,InsertMode mode,Vec l) 32 { 33 PetscErrorCode ierr; 34 DM_Shell *shell = (DM_Shell*)dm->data; 35 36 PetscFunctionBegin; 37 if (!shell->gtol) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE, "Cannot be used without first setting the scatter context via DMShellSetGlobalToLocalVecScatter()"); 38 ierr = VecScatterBegin(shell->gtol,g,l,mode,SCATTER_FORWARD);CHKERRQ(ierr); 39 PetscFunctionReturn(0); 40 } 41 42 /*@ 43 DMGlobalToLocalEndDefaultShell - Uses the GlobalToLocal VecScatter context set by the user to end a global to local scatter 44 Collective 45 46 Input Arguments: 47 + dm - shell DM 48 . g - global vector 49 . mode - InsertMode 50 - l - local vector 51 52 Level: advanced 53 54 .seealso: DMGlobalToLocalBeginDefaultShell() 55 @*/ 56 PetscErrorCode DMGlobalToLocalEndDefaultShell(DM dm,Vec g,InsertMode mode,Vec l) 57 { 58 PetscErrorCode ierr; 59 DM_Shell *shell = (DM_Shell*)dm->data; 60 61 PetscFunctionBegin; 62 if (!shell->gtol) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE, "Cannot be used without first setting the scatter context via DMShellSetGlobalToLocalVecScatter()"); 63 ierr = VecScatterEnd(shell->gtol,g,l,mode,SCATTER_FORWARD);CHKERRQ(ierr); 64 PetscFunctionReturn(0); 65 } 66 67 /*@ 68 DMLocalToGlobalBeginDefaultShell - Uses the LocalToGlobal VecScatter context set by the user to begin a local to global scatter 69 Collective 70 71 Input Arguments: 72 + dm - shell DM 73 . l - local vector 74 . mode - InsertMode 75 - g - global vector 76 77 Level: advanced 78 79 Note: This is not normally called directly by user code, generally user code calls DMLocalToGlobalBegin() and DMLocalToGlobalEnd(). If the user provides their own custom routines to DMShellSetLocalToGlobal() then those routines might have reason to call this function. 80 81 .seealso: DMLocalToGlobalEndDefaultShell() 82 @*/ 83 PetscErrorCode DMLocalToGlobalBeginDefaultShell(DM dm,Vec l,InsertMode mode,Vec g) 84 { 85 PetscErrorCode ierr; 86 DM_Shell *shell = (DM_Shell*)dm->data; 87 88 PetscFunctionBegin; 89 if (!shell->ltog) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE, "Cannot be used without first setting the scatter context via DMShellSetLocalToGlobalVecScatter()"); 90 ierr = VecScatterBegin(shell->ltog,l,g,mode,SCATTER_FORWARD);CHKERRQ(ierr); 91 PetscFunctionReturn(0); 92 } 93 94 /*@ 95 DMLocalToGlobalEndDefaultShell - Uses the LocalToGlobal VecScatter context set by the user to end a local to global scatter 96 Collective 97 98 Input Arguments: 99 + dm - shell DM 100 . l - local vector 101 . mode - InsertMode 102 - g - global vector 103 104 Level: advanced 105 106 .seealso: DMLocalToGlobalBeginDefaultShell() 107 @*/ 108 PetscErrorCode DMLocalToGlobalEndDefaultShell(DM dm,Vec l,InsertMode mode,Vec g) 109 { 110 PetscErrorCode ierr; 111 DM_Shell *shell = (DM_Shell*)dm->data; 112 113 PetscFunctionBegin; 114 if (!shell->ltog) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE, "Cannot be used without first setting the scatter context via DMShellSetLocalToGlobalVecScatter()"); 115 ierr = VecScatterEnd(shell->ltog,l,g,mode,SCATTER_FORWARD);CHKERRQ(ierr); 116 PetscFunctionReturn(0); 117 } 118 119 /*@ 120 DMLocalToLocalBeginDefaultShell - Uses the LocalToLocal VecScatter context set by the user to begin a local to local scatter 121 Collective 122 123 Input Arguments: 124 + dm - shell DM 125 . g - the original local vector 126 - mode - InsertMode 127 128 Output Parameter: 129 . l - the local vector with correct ghost values 130 131 Level: advanced 132 133 Note: This is not normally called directly by user code, generally user code calls DMLocalToLocalBegin() and DMLocalToLocalEnd(). If the user provides their own custom routines to DMShellSetLocalToLocal() then those routines might have reason to call this function. 134 135 .seealso: DMLocalToLocalEndDefaultShell() 136 @*/ 137 PetscErrorCode DMLocalToLocalBeginDefaultShell(DM dm,Vec g,InsertMode mode,Vec l) 138 { 139 PetscErrorCode ierr; 140 DM_Shell *shell = (DM_Shell*)dm->data; 141 142 PetscFunctionBegin; 143 if (!shell->ltol) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE, "Cannot be used without first setting the scatter context via DMShellSetLocalToLocalVecScatter()"); 144 ierr = VecScatterBegin(shell->ltol,g,l,mode,SCATTER_FORWARD);CHKERRQ(ierr); 145 PetscFunctionReturn(0); 146 } 147 148 /*@ 149 DMLocalToLocalEndDefaultShell - Uses the LocalToLocal VecScatter context set by the user to end a local to local scatter 150 Collective 151 152 Input Arguments: 153 + dm - shell DM 154 . g - the original local vector 155 - mode - InsertMode 156 157 Output Parameter: 158 . l - the local vector with correct ghost values 159 160 Level: advanced 161 162 .seealso: DMLocalToLocalBeginDefaultShell() 163 @*/ 164 PetscErrorCode DMLocalToLocalEndDefaultShell(DM dm,Vec g,InsertMode mode,Vec l) 165 { 166 PetscErrorCode ierr; 167 DM_Shell *shell = (DM_Shell*)dm->data; 168 169 PetscFunctionBegin; 170 if (!shell->ltol) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE, "Cannot be used without first setting the scatter context via DMShellSetGlobalToLocalVecScatter()"); 171 ierr = VecScatterEnd(shell->ltol,g,l,mode,SCATTER_FORWARD);CHKERRQ(ierr); 172 PetscFunctionReturn(0); 173 } 174 175 static PetscErrorCode DMCreateMatrix_Shell(DM dm,Mat *J) 176 { 177 PetscErrorCode ierr; 178 DM_Shell *shell = (DM_Shell*)dm->data; 179 Mat A; 180 181 PetscFunctionBegin; 182 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 183 PetscValidPointer(J,3); 184 if (!shell->A) { 185 if (shell->Xglobal) { 186 PetscInt m,M; 187 ierr = PetscInfo(dm,"Naively creating matrix using global vector distribution without preallocation\n");CHKERRQ(ierr); 188 ierr = VecGetSize(shell->Xglobal,&M);CHKERRQ(ierr); 189 ierr = VecGetLocalSize(shell->Xglobal,&m);CHKERRQ(ierr); 190 ierr = MatCreate(PetscObjectComm((PetscObject)dm),&shell->A);CHKERRQ(ierr); 191 ierr = MatSetSizes(shell->A,m,m,M,M);CHKERRQ(ierr); 192 ierr = MatSetType(shell->A,dm->mattype);CHKERRQ(ierr); 193 ierr = MatSetUp(shell->A);CHKERRQ(ierr); 194 } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Must call DMShellSetMatrix(), DMShellSetCreateMatrix(), or provide a vector"); 195 } 196 A = shell->A; 197 /* the check below is tacky and incomplete */ 198 if (dm->mattype) { 199 PetscBool flg,aij,seqaij,mpiaij; 200 ierr = PetscObjectTypeCompare((PetscObject)A,dm->mattype,&flg);CHKERRQ(ierr); 201 ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&seqaij);CHKERRQ(ierr); 202 ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&mpiaij);CHKERRQ(ierr); 203 ierr = PetscStrcmp(dm->mattype,MATAIJ,&aij);CHKERRQ(ierr); 204 if (!flg) { 205 if (!(aij && (seqaij || mpiaij))) SETERRQ2(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_NOTSAMETYPE,"Requested matrix of type %s, but only %s available",dm->mattype,((PetscObject)A)->type_name); 206 } 207 } 208 if (((PetscObject)A)->refct < 2) { /* We have an exclusive reference so we can give it out */ 209 ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 210 ierr = MatZeroEntries(A);CHKERRQ(ierr); 211 *J = A; 212 } else { /* Need to create a copy, could use MAT_SHARE_NONZERO_PATTERN in most cases */ 213 ierr = MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,J);CHKERRQ(ierr); 214 ierr = MatZeroEntries(*J);CHKERRQ(ierr); 215 } 216 PetscFunctionReturn(0); 217 } 218 219 PetscErrorCode DMCreateGlobalVector_Shell(DM dm,Vec *gvec) 220 { 221 PetscErrorCode ierr; 222 DM_Shell *shell = (DM_Shell*)dm->data; 223 Vec X; 224 225 PetscFunctionBegin; 226 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 227 PetscValidPointer(gvec,2); 228 *gvec = 0; 229 X = shell->Xglobal; 230 if (!X) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Must call DMShellSetGlobalVector() or DMShellSetCreateGlobalVector()"); 231 if (((PetscObject)X)->refct < 2) { /* We have an exclusive reference so we can give it out */ 232 ierr = PetscObjectReference((PetscObject)X);CHKERRQ(ierr); 233 ierr = VecZeroEntries(X);CHKERRQ(ierr); 234 *gvec = X; 235 } else { /* Need to create a copy, could use MAT_SHARE_NONZERO_PATTERN in most cases */ 236 ierr = VecDuplicate(X,gvec);CHKERRQ(ierr); 237 ierr = VecZeroEntries(*gvec);CHKERRQ(ierr); 238 } 239 ierr = VecSetDM(*gvec,dm);CHKERRQ(ierr); 240 PetscFunctionReturn(0); 241 } 242 243 PetscErrorCode DMCreateLocalVector_Shell(DM dm,Vec *gvec) 244 { 245 PetscErrorCode ierr; 246 DM_Shell *shell = (DM_Shell*)dm->data; 247 Vec X; 248 249 PetscFunctionBegin; 250 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 251 PetscValidPointer(gvec,2); 252 *gvec = 0; 253 X = shell->Xlocal; 254 if (!X) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Must call DMShellSetLocalVector() or DMShellSetCreateLocalVector()"); 255 if (((PetscObject)X)->refct < 2) { /* We have an exclusive reference so we can give it out */ 256 ierr = PetscObjectReference((PetscObject)X);CHKERRQ(ierr); 257 ierr = VecZeroEntries(X);CHKERRQ(ierr); 258 *gvec = X; 259 } else { /* Need to create a copy, could use MAT_SHARE_NONZERO_PATTERN in most cases */ 260 ierr = VecDuplicate(X,gvec);CHKERRQ(ierr); 261 ierr = VecZeroEntries(*gvec);CHKERRQ(ierr); 262 } 263 ierr = VecSetDM(*gvec,dm);CHKERRQ(ierr); 264 PetscFunctionReturn(0); 265 } 266 267 /*@ 268 DMShellSetContext - set some data to be usable by this DM 269 270 Collective 271 272 Input Arguments: 273 + dm - shell DM 274 - ctx - the context 275 276 Level: advanced 277 278 .seealso: DMCreateMatrix(), DMShellGetContext() 279 @*/ 280 PetscErrorCode DMShellSetContext(DM dm,void *ctx) 281 { 282 DM_Shell *shell = (DM_Shell*)dm->data; 283 PetscErrorCode ierr; 284 PetscBool isshell; 285 286 PetscFunctionBegin; 287 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 288 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 289 if (!isshell) PetscFunctionReturn(0); 290 shell->ctx = ctx; 291 PetscFunctionReturn(0); 292 } 293 294 /*@ 295 DMShellGetContext - set some data to be usable by this DM 296 297 Collective 298 299 Input Argument: 300 . dm - shell DM 301 302 Output Argument: 303 . ctx - the context 304 305 Level: advanced 306 307 .seealso: DMCreateMatrix(), DMShellSetContext() 308 @*/ 309 PetscErrorCode DMShellGetContext(DM dm,void **ctx) 310 { 311 DM_Shell *shell = (DM_Shell*)dm->data; 312 PetscErrorCode ierr; 313 PetscBool isshell; 314 315 PetscFunctionBegin; 316 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 317 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 318 if (!isshell) PetscFunctionReturn(0); 319 *ctx = shell->ctx; 320 PetscFunctionReturn(0); 321 } 322 323 /*@ 324 DMShellSetMatrix - sets a template matrix associated with the DMShell 325 326 Collective 327 328 Input Arguments: 329 + dm - shell DM 330 - J - template matrix 331 332 Level: advanced 333 334 .seealso: DMCreateMatrix(), DMShellSetCreateMatrix(), DMShellSetContext(), DMShellGetContext() 335 @*/ 336 PetscErrorCode DMShellSetMatrix(DM dm,Mat J) 337 { 338 DM_Shell *shell = (DM_Shell*)dm->data; 339 PetscErrorCode ierr; 340 PetscBool isshell; 341 342 PetscFunctionBegin; 343 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 344 PetscValidHeaderSpecific(J,MAT_CLASSID,2); 345 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 346 if (!isshell) PetscFunctionReturn(0); 347 ierr = PetscObjectReference((PetscObject)J);CHKERRQ(ierr); 348 ierr = MatDestroy(&shell->A);CHKERRQ(ierr); 349 shell->A = J; 350 PetscFunctionReturn(0); 351 } 352 353 /*@C 354 DMShellSetCreateMatrix - sets the routine to create a matrix associated with the shell DM 355 356 Logically Collective on DM 357 358 Input Arguments: 359 + dm - the shell DM 360 - func - the function to create a matrix 361 362 Level: advanced 363 364 .seealso: DMCreateMatrix(), DMShellSetMatrix(), DMShellSetContext(), DMShellGetContext() 365 @*/ 366 PetscErrorCode DMShellSetCreateMatrix(DM dm,PetscErrorCode (*func)(DM,Mat*)) 367 { 368 369 PetscFunctionBegin; 370 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 371 dm->ops->creatematrix = func; 372 PetscFunctionReturn(0); 373 } 374 375 /*@ 376 DMShellSetGlobalVector - sets a template global vector associated with the DMShell 377 378 Logically Collective on DM 379 380 Input Arguments: 381 + dm - shell DM 382 - X - template vector 383 384 Level: advanced 385 386 .seealso: DMCreateGlobalVector(), DMShellSetMatrix(), DMShellSetCreateGlobalVector() 387 @*/ 388 PetscErrorCode DMShellSetGlobalVector(DM dm,Vec X) 389 { 390 DM_Shell *shell = (DM_Shell*)dm->data; 391 PetscErrorCode ierr; 392 PetscBool isshell; 393 DM vdm; 394 395 PetscFunctionBegin; 396 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 397 PetscValidHeaderSpecific(X,VEC_CLASSID,2); 398 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 399 if (!isshell) PetscFunctionReturn(0); 400 ierr = VecGetDM(X,&vdm);CHKERRQ(ierr); 401 /* 402 if the vector proposed as the new base global vector for the DM is a DM vector associated 403 with the same DM then the current base global vector for the DM is ok and if we replace it with the new one 404 we get a circular dependency that prevents the DM from being destroy when it should be. 405 This occurs when SNESSet/GetNPC() is used with a SNES that does not have a user provided 406 DM attached to it since the inner SNES (which shares the DM with the outer SNES) tries 407 to set its input vector (which is associated with the DM) as the base global vector. 408 Thanks to Juan P. Mendez Granado Re: [petsc-maint] Nonlinear conjugate gradien 409 for pointing out the problem. 410 */ 411 if (vdm == dm) PetscFunctionReturn(0); 412 ierr = PetscObjectReference((PetscObject)X);CHKERRQ(ierr); 413 ierr = VecDestroy(&shell->Xglobal);CHKERRQ(ierr); 414 shell->Xglobal = X; 415 PetscFunctionReturn(0); 416 } 417 418 /*@C 419 DMShellSetCreateGlobalVector - sets the routine to create a global vector associated with the shell DM 420 421 Logically Collective 422 423 Input Arguments: 424 + dm - the shell DM 425 - func - the creation routine 426 427 Level: advanced 428 429 .seealso: DMShellSetGlobalVector(), DMShellSetCreateMatrix(), DMShellSetContext(), DMShellGetContext() 430 @*/ 431 PetscErrorCode DMShellSetCreateGlobalVector(DM dm,PetscErrorCode (*func)(DM,Vec*)) 432 { 433 434 PetscFunctionBegin; 435 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 436 dm->ops->createglobalvector = func; 437 PetscFunctionReturn(0); 438 } 439 440 /*@ 441 DMShellSetLocalVector - sets a template local vector associated with the DMShell 442 443 Logically Collective on DM 444 445 Input Arguments: 446 + dm - shell DM 447 - X - template vector 448 449 Level: advanced 450 451 .seealso: DMCreateLocalVector(), DMShellSetMatrix(), DMShellSetCreateLocalVector() 452 @*/ 453 PetscErrorCode DMShellSetLocalVector(DM dm,Vec X) 454 { 455 DM_Shell *shell = (DM_Shell*)dm->data; 456 PetscErrorCode ierr; 457 PetscBool isshell; 458 DM vdm; 459 460 PetscFunctionBegin; 461 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 462 PetscValidHeaderSpecific(X,VEC_CLASSID,2); 463 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 464 if (!isshell) PetscFunctionReturn(0); 465 ierr = VecGetDM(X,&vdm);CHKERRQ(ierr); 466 /* 467 if the vector proposed as the new base global vector for the DM is a DM vector associated 468 with the same DM then the current base global vector for the DM is ok and if we replace it with the new one 469 we get a circular dependency that prevents the DM from being destroy when it should be. 470 This occurs when SNESSet/GetNPC() is used with a SNES that does not have a user provided 471 DM attached to it since the inner SNES (which shares the DM with the outer SNES) tries 472 to set its input vector (which is associated with the DM) as the base global vector. 473 Thanks to Juan P. Mendez Granado Re: [petsc-maint] Nonlinear conjugate gradien 474 for pointing out the problem. 475 */ 476 if (vdm == dm) PetscFunctionReturn(0); 477 ierr = PetscObjectReference((PetscObject)X);CHKERRQ(ierr); 478 ierr = VecDestroy(&shell->Xlocal);CHKERRQ(ierr); 479 shell->Xlocal = X; 480 PetscFunctionReturn(0); 481 } 482 483 /*@C 484 DMShellSetCreateLocalVector - sets the routine to create a local vector associated with the shell DM 485 486 Logically Collective 487 488 Input Arguments: 489 + dm - the shell DM 490 - func - the creation routine 491 492 Level: advanced 493 494 .seealso: DMShellSetLocalVector(), DMShellSetCreateMatrix(), DMShellSetContext(), DMShellGetContext() 495 @*/ 496 PetscErrorCode DMShellSetCreateLocalVector(DM dm,PetscErrorCode (*func)(DM,Vec*)) 497 { 498 499 PetscFunctionBegin; 500 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 501 dm->ops->createlocalvector = func; 502 PetscFunctionReturn(0); 503 } 504 505 /*@C 506 DMShellSetGlobalToLocal - Sets the routines used to perform a global to local scatter 507 508 Logically Collective on DM 509 510 Input Arguments 511 + dm - the shell DM 512 . begin - the routine that begins the global to local scatter 513 - end - the routine that ends the global to local scatter 514 515 Notes: 516 If these functions are not provided but DMShellSetGlobalToLocalVecScatter() is called then 517 DMGlobalToLocalBeginDefaultShell()/DMGlobalToLocalEndDefaultShell() are used to to perform the transfers 518 519 Level: advanced 520 521 .seealso: DMShellSetLocalToGlobal(), DMGlobalToLocalBeginDefaultShell(), DMGlobalToLocalEndDefaultShell() 522 @*/ 523 PetscErrorCode DMShellSetGlobalToLocal(DM dm,PetscErrorCode (*begin)(DM,Vec,InsertMode,Vec),PetscErrorCode (*end)(DM,Vec,InsertMode,Vec)) { 524 PetscFunctionBegin; 525 dm->ops->globaltolocalbegin = begin; 526 dm->ops->globaltolocalend = end; 527 PetscFunctionReturn(0); 528 } 529 530 /*@C 531 DMShellSetLocalToGlobal - Sets the routines used to perform a local to global scatter 532 533 Logically Collective on DM 534 535 Input Arguments 536 + dm - the shell DM 537 . begin - the routine that begins the local to global scatter 538 - end - the routine that ends the local to global scatter 539 540 Notes: 541 If these functions are not provided but DMShellSetLocalToGlobalVecScatter() is called then 542 DMLocalToGlobalBeginDefaultShell()/DMLocalToGlobalEndDefaultShell() are used to to perform the transfers 543 544 Level: advanced 545 546 .seealso: DMShellSetGlobalToLocal() 547 @*/ 548 PetscErrorCode DMShellSetLocalToGlobal(DM dm,PetscErrorCode (*begin)(DM,Vec,InsertMode,Vec),PetscErrorCode (*end)(DM,Vec,InsertMode,Vec)) { 549 PetscFunctionBegin; 550 dm->ops->localtoglobalbegin = begin; 551 dm->ops->localtoglobalend = end; 552 PetscFunctionReturn(0); 553 } 554 555 /*@C 556 DMShellSetLocalToLocal - Sets the routines used to perform a local to local scatter 557 558 Logically Collective on DM 559 560 Input Arguments 561 + dm - the shell DM 562 . begin - the routine that begins the local to local scatter 563 - end - the routine that ends the local to local scatter 564 565 Notes: 566 If these functions are not provided but DMShellSetLocalToLocalVecScatter() is called then 567 DMLocalToLocalBeginDefaultShell()/DMLocalToLocalEndDefaultShell() are used to to perform the transfers 568 569 Level: advanced 570 571 .seealso: DMShellSetGlobalToLocal(), DMLocalToLocalBeginDefaultShell(), DMLocalToLocalEndDefaultShell() 572 @*/ 573 PetscErrorCode DMShellSetLocalToLocal(DM dm,PetscErrorCode (*begin)(DM,Vec,InsertMode,Vec),PetscErrorCode (*end)(DM,Vec,InsertMode,Vec)) { 574 PetscFunctionBegin; 575 dm->ops->localtolocalbegin = begin; 576 dm->ops->localtolocalend = end; 577 PetscFunctionReturn(0); 578 } 579 580 /*@ 581 DMShellSetGlobalToLocalVecScatter - Sets a VecScatter context for global to local communication 582 583 Logically Collective on DM 584 585 Input Arguments 586 + dm - the shell DM 587 - gtol - the global to local VecScatter context 588 589 Level: advanced 590 591 .seealso: DMShellSetGlobalToLocal(), DMGlobalToLocalBeginDefaultShell(), DMGlobalToLocalEndDefaultShell() 592 @*/ 593 PetscErrorCode DMShellSetGlobalToLocalVecScatter(DM dm, VecScatter gtol) 594 { 595 DM_Shell *shell = (DM_Shell*)dm->data; 596 PetscErrorCode ierr; 597 598 PetscFunctionBegin; 599 ierr = PetscObjectReference((PetscObject)gtol);CHKERRQ(ierr); 600 /* Call VecScatterDestroy() to avoid a memory leak in case of re-setting. */ 601 ierr = VecScatterDestroy(&shell->gtol);CHKERRQ(ierr); 602 shell->gtol = gtol; 603 PetscFunctionReturn(0); 604 } 605 606 /*@ 607 DMShellSetLocalToGlobalVecScatter - Sets a VecScatter context for local to global communication 608 609 Logically Collective on DM 610 611 Input Arguments 612 + dm - the shell DM 613 - ltog - the local to global VecScatter context 614 615 Level: advanced 616 617 .seealso: DMShellSetLocalToGlobal(), DMLocalToGlobalBeginDefaultShell(), DMLocalToGlobalEndDefaultShell() 618 @*/ 619 PetscErrorCode DMShellSetLocalToGlobalVecScatter(DM dm, VecScatter ltog) 620 { 621 DM_Shell *shell = (DM_Shell*)dm->data; 622 PetscErrorCode ierr; 623 624 PetscFunctionBegin; 625 ierr = PetscObjectReference((PetscObject)ltog);CHKERRQ(ierr); 626 /* Call VecScatterDestroy() to avoid a memory leak in case of re-setting. */ 627 ierr = VecScatterDestroy(&shell->ltog);CHKERRQ(ierr); 628 shell->ltog = ltog; 629 PetscFunctionReturn(0); 630 } 631 632 /*@ 633 DMShellSetLocalToLocalVecScatter - Sets a VecScatter context for local to local communication 634 635 Logically Collective on DM 636 637 Input Arguments 638 + dm - the shell DM 639 - ltol - the local to local VecScatter context 640 641 Level: advanced 642 643 .seealso: DMShellSetLocalToLocal(), DMLocalToLocalBeginDefaultShell(), DMLocalToLocalEndDefaultShell() 644 @*/ 645 PetscErrorCode DMShellSetLocalToLocalVecScatter(DM dm, VecScatter ltol) 646 { 647 DM_Shell *shell = (DM_Shell*)dm->data; 648 PetscErrorCode ierr; 649 650 PetscFunctionBegin; 651 ierr = PetscObjectReference((PetscObject)ltol);CHKERRQ(ierr); 652 /* Call VecScatterDestroy() to avoid a memory leak in case of re-setting. */ 653 ierr = VecScatterDestroy(&shell->ltol);CHKERRQ(ierr); 654 shell->ltol = ltol; 655 PetscFunctionReturn(0); 656 } 657 658 /*@C 659 DMShellSetCoarsen - Set the routine used to coarsen the shell DM 660 661 Logically Collective on DM 662 663 Input Arguments 664 + dm - the shell DM 665 - coarsen - the routine that coarsens the DM 666 667 Level: advanced 668 669 .seealso: DMShellSetRefine(), DMCoarsen(), DMShellSetContext(), DMShellGetContext() 670 @*/ 671 PetscErrorCode DMShellSetCoarsen(DM dm, PetscErrorCode (*coarsen)(DM,MPI_Comm,DM*)) 672 { 673 PetscErrorCode ierr; 674 PetscBool isshell; 675 676 PetscFunctionBegin; 677 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 678 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 679 if (!isshell) PetscFunctionReturn(0); 680 dm->ops->coarsen = coarsen; 681 PetscFunctionReturn(0); 682 } 683 684 /*@C 685 DMShellSetRefine - Set the routine used to refine the shell DM 686 687 Logically Collective on DM 688 689 Input Arguments 690 + dm - the shell DM 691 - refine - the routine that refines the DM 692 693 Level: advanced 694 695 .seealso: DMShellSetCoarsen(), DMRefine(), DMShellSetContext(), DMShellGetContext() 696 @*/ 697 PetscErrorCode DMShellSetRefine(DM dm, PetscErrorCode (*refine)(DM,MPI_Comm,DM*)) 698 { 699 PetscErrorCode ierr; 700 PetscBool isshell; 701 702 PetscFunctionBegin; 703 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 704 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 705 if (!isshell) PetscFunctionReturn(0); 706 dm->ops->refine = refine; 707 PetscFunctionReturn(0); 708 } 709 710 /*@C 711 DMShellSetCreateInterpolation - Set the routine used to create the interpolation operator 712 713 Logically Collective on DM 714 715 Input Arguments 716 + dm - the shell DM 717 - interp - the routine to create the interpolation 718 719 Level: advanced 720 721 .seealso: DMShellSetCreateInjection(), DMCreateInterpolation(), DMShellSetCreateRestriction(), DMShellSetContext(), DMShellGetContext() 722 @*/ 723 PetscErrorCode DMShellSetCreateInterpolation(DM dm, PetscErrorCode (*interp)(DM,DM,Mat*,Vec*)) 724 { 725 PetscErrorCode ierr; 726 PetscBool isshell; 727 728 PetscFunctionBegin; 729 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 730 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 731 if (!isshell) PetscFunctionReturn(0); 732 dm->ops->createinterpolation = interp; 733 PetscFunctionReturn(0); 734 } 735 736 /*@C 737 DMShellSetCreateRestriction - Set the routine used to create the restriction operator 738 739 Logically Collective on DM 740 741 Input Arguments 742 + dm - the shell DM 743 - striction- the routine to create the restriction 744 745 Level: advanced 746 747 .seealso: DMShellSetCreateInjection(), DMCreateInterpolation(), DMShellSetContext(), DMShellGetContext() 748 @*/ 749 PetscErrorCode DMShellSetCreateRestriction(DM dm, PetscErrorCode (*restriction)(DM,DM,Mat*)) 750 { 751 PetscErrorCode ierr; 752 PetscBool isshell; 753 754 PetscFunctionBegin; 755 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 756 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 757 if (!isshell) PetscFunctionReturn(0); 758 dm->ops->createrestriction = restriction; 759 PetscFunctionReturn(0); 760 } 761 762 /*@C 763 DMShellSetCreateInjection - Set the routine used to create the injection operator 764 765 Logically Collective on DM 766 767 Input Arguments 768 + dm - the shell DM 769 - inject - the routine to create the injection 770 771 Level: advanced 772 773 .seealso: DMShellSetCreateInterpolation(), DMCreateInjection(), DMShellSetContext(), DMShellGetContext() 774 @*/ 775 PetscErrorCode DMShellSetCreateInjection(DM dm, PetscErrorCode (*inject)(DM,DM,Mat*)) 776 { 777 PetscErrorCode ierr; 778 PetscBool isshell; 779 780 PetscFunctionBegin; 781 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 782 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 783 if (!isshell) PetscFunctionReturn(0); 784 dm->ops->getinjection = inject; 785 PetscFunctionReturn(0); 786 } 787 788 static PetscErrorCode DMHasCreateInjection_Shell(DM dm, PetscBool *flg) 789 { 790 PetscErrorCode ierr; 791 PetscBool isshell; 792 793 PetscFunctionBegin; 794 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 795 PetscValidPointer(flg,2); 796 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 797 if (!isshell) { 798 *flg = PETSC_FALSE; 799 } else { 800 *flg = dm->ops->getinjection ? PETSC_TRUE : PETSC_FALSE; 801 } 802 PetscFunctionReturn(0); 803 } 804 /*@C 805 DMShellSetCreateFieldDecomposition - Set the routine used to create a decomposition of fields for the shell DM 806 807 Logically Collective on DM 808 809 Input Arguments 810 + dm - the shell DM 811 - decomp - the routine to create the decomposition 812 813 Level: advanced 814 815 .seealso: DMCreateFieldDecomposition(), DMShellSetContext(), DMShellGetContext() 816 @*/ 817 PetscErrorCode DMShellSetCreateFieldDecomposition(DM dm, PetscErrorCode (*decomp)(DM,PetscInt*,char***, IS**,DM**)) 818 { 819 PetscErrorCode ierr; 820 PetscBool isshell; 821 822 PetscFunctionBegin; 823 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 824 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 825 if (!isshell) PetscFunctionReturn(0); 826 dm->ops->createfielddecomposition = decomp; 827 PetscFunctionReturn(0); 828 } 829 830 /*@C 831 DMShellSetCreateDomainDecomposition - Set the routine used to create a domain decomposition for the shell DM 832 833 Logically Collective on DM 834 835 Input Arguments 836 + dm - the shell DM 837 - decomp - the routine to create the decomposition 838 839 Level: advanced 840 841 .seealso: DMCreateDomainDecomposition(), DMShellSetContext(), DMShellGetContext() 842 @*/ 843 PetscErrorCode DMShellSetCreateDomainDecomposition(DM dm, PetscErrorCode (*decomp)(DM,PetscInt*,char***, IS**,IS**,DM**)) 844 { 845 PetscErrorCode ierr; 846 PetscBool isshell; 847 848 PetscFunctionBegin; 849 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 850 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 851 if (!isshell) PetscFunctionReturn(0); 852 dm->ops->createdomaindecomposition = decomp; 853 PetscFunctionReturn(0); 854 } 855 856 /*@C 857 DMShellSetCreateDomainDecompositionScatters - Set the routine used to create the scatter contexts for domain decomposition with a shell DM 858 859 Logically Collective on DM 860 861 Input Arguments 862 + dm - the shell DM 863 - scatter - the routine to create the scatters 864 865 Level: advanced 866 867 .seealso: DMCreateDomainDecompositionScatters(), DMShellSetContext(), DMShellGetContext() 868 @*/ 869 PetscErrorCode DMShellSetCreateDomainDecompositionScatters(DM dm, PetscErrorCode (*scatter)(DM,PetscInt,DM*,VecScatter**,VecScatter**,VecScatter**)) 870 { 871 PetscErrorCode ierr; 872 PetscBool isshell; 873 874 PetscFunctionBegin; 875 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 876 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 877 if (!isshell) PetscFunctionReturn(0); 878 dm->ops->createddscatters = scatter; 879 PetscFunctionReturn(0); 880 } 881 882 /*@C 883 DMShellSetCreateSubDM - Set the routine used to create a sub DM from the shell DM 884 885 Logically Collective on DM 886 887 Input Arguments 888 + dm - the shell DM 889 - subdm - the routine to create the decomposition 890 891 Level: advanced 892 893 .seealso: DMCreateSubDM(), DMShellSetContext(), DMShellGetContext() 894 @*/ 895 PetscErrorCode DMShellSetCreateSubDM(DM dm, PetscErrorCode (*subdm)(DM,PetscInt,const PetscInt[],IS*,DM*)) 896 { 897 PetscErrorCode ierr; 898 PetscBool isshell; 899 900 PetscFunctionBegin; 901 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 902 ierr = PetscObjectTypeCompare((PetscObject)dm,DMSHELL,&isshell);CHKERRQ(ierr); 903 if (!isshell) PetscFunctionReturn(0); 904 dm->ops->createsubdm = subdm; 905 PetscFunctionReturn(0); 906 } 907 908 static PetscErrorCode DMDestroy_Shell(DM dm) 909 { 910 PetscErrorCode ierr; 911 DM_Shell *shell = (DM_Shell*)dm->data; 912 913 PetscFunctionBegin; 914 ierr = MatDestroy(&shell->A);CHKERRQ(ierr); 915 ierr = VecDestroy(&shell->Xglobal);CHKERRQ(ierr); 916 ierr = VecDestroy(&shell->Xlocal);CHKERRQ(ierr); 917 ierr = VecScatterDestroy(&shell->gtol);CHKERRQ(ierr); 918 ierr = VecScatterDestroy(&shell->ltog);CHKERRQ(ierr); 919 ierr = VecScatterDestroy(&shell->ltol);CHKERRQ(ierr); 920 /* This was originally freed in DMDestroy(), but that prevents reference counting of backend objects */ 921 ierr = PetscFree(shell);CHKERRQ(ierr); 922 PetscFunctionReturn(0); 923 } 924 925 static PetscErrorCode DMView_Shell(DM dm,PetscViewer v) 926 { 927 PetscErrorCode ierr; 928 DM_Shell *shell = (DM_Shell*)dm->data; 929 930 PetscFunctionBegin; 931 ierr = VecView(shell->Xglobal,v);CHKERRQ(ierr); 932 PetscFunctionReturn(0); 933 } 934 935 static PetscErrorCode DMLoad_Shell(DM dm,PetscViewer v) 936 { 937 PetscErrorCode ierr; 938 DM_Shell *shell = (DM_Shell*)dm->data; 939 940 PetscFunctionBegin; 941 ierr = VecCreate(PetscObjectComm((PetscObject)dm),&shell->Xglobal);CHKERRQ(ierr); 942 ierr = VecLoad(shell->Xglobal,v);CHKERRQ(ierr); 943 PetscFunctionReturn(0); 944 } 945 946 PetscErrorCode DMCreateSubDM_Shell(DM dm, PetscInt numFields, const PetscInt fields[], IS *is, DM *subdm) 947 { 948 PetscErrorCode ierr; 949 950 PetscFunctionBegin; 951 if (subdm) {ierr = DMShellCreate(PetscObjectComm((PetscObject) dm), subdm);CHKERRQ(ierr);} 952 ierr = DMCreateSubDM_Section_Private(dm, numFields, fields, is, subdm);CHKERRQ(ierr); 953 PetscFunctionReturn(0); 954 } 955 956 PETSC_EXTERN PetscErrorCode DMCreate_Shell(DM dm) 957 { 958 PetscErrorCode ierr; 959 DM_Shell *shell; 960 961 PetscFunctionBegin; 962 ierr = PetscNewLog(dm,&shell);CHKERRQ(ierr); 963 dm->data = shell; 964 965 ierr = PetscObjectChangeTypeName((PetscObject)dm,DMSHELL);CHKERRQ(ierr); 966 967 dm->ops->destroy = DMDestroy_Shell; 968 dm->ops->createglobalvector = DMCreateGlobalVector_Shell; 969 dm->ops->createlocalvector = DMCreateLocalVector_Shell; 970 dm->ops->creatematrix = DMCreateMatrix_Shell; 971 dm->ops->view = DMView_Shell; 972 dm->ops->load = DMLoad_Shell; 973 dm->ops->globaltolocalbegin = DMGlobalToLocalBeginDefaultShell; 974 dm->ops->globaltolocalend = DMGlobalToLocalEndDefaultShell; 975 dm->ops->localtoglobalbegin = DMLocalToGlobalBeginDefaultShell; 976 dm->ops->localtoglobalend = DMLocalToGlobalEndDefaultShell; 977 dm->ops->localtolocalbegin = DMLocalToLocalBeginDefaultShell; 978 dm->ops->localtolocalend = DMLocalToLocalEndDefaultShell; 979 dm->ops->createsubdm = DMCreateSubDM_Shell; 980 dm->ops->hascreateinjection = DMHasCreateInjection_Shell; 981 PetscFunctionReturn(0); 982 } 983 984 /*@ 985 DMShellCreate - Creates a shell DM object, used to manage user-defined problem data 986 987 Collective on MPI_Comm 988 989 Input Parameter: 990 . comm - the processors that will share the global vector 991 992 Output Parameters: 993 . shell - the shell DM 994 995 Level: advanced 996 997 .seealso DMDestroy(), DMCreateGlobalVector(), DMCreateLocalVector(), DMShellSetContext(), DMShellGetContext() 998 @*/ 999 PetscErrorCode DMShellCreate(MPI_Comm comm,DM *dm) 1000 { 1001 PetscErrorCode ierr; 1002 1003 PetscFunctionBegin; 1004 PetscValidPointer(dm,2); 1005 ierr = DMCreate(comm,dm);CHKERRQ(ierr); 1006 ierr = DMSetType(*dm,DMSHELL);CHKERRQ(ierr); 1007 ierr = DMSetUp(*dm);CHKERRQ(ierr); 1008 PetscFunctionReturn(0); 1009 } 1010 1011