1 /* TODOLIST 2 3 ConstraintsSetup 4 - tolerances for constraints as an option (take care of single precision!) 5 - Can MAT_IGNORE_ZERO_ENTRIES be used for Constraints Matrix? 6 7 Solvers 8 - Add support for reuse fill and cholecky factor for coarse solver (similar to local solvers) 9 - Propagate ksp prefixes for solvers to mat objects? 10 - Propagate nearnullspace info among levels 11 12 User interface 13 - Negative indices in dirichlet and Neumann is should be skipped (now they cause out-of-bounds access) 14 - Provide PCApplyTranpose_BDDC 15 - DofSplitting and DM attached to pc? 16 17 Debugging output 18 - Better management of verbosity levels of debugging output 19 20 Build 21 - make runexe59 22 23 Extra 24 - Is it possible to work with PCBDDCGraph on boundary indices only (less memory consumed)? 25 - Why options for "pc_bddc_coarse" solver gets propagated to "pc_bddc_coarse_1" solver? 26 - add support for computing h,H and related using coordinates? 27 - Change of basis approach does not work with my nonlinear mechanics example. why? (seems not an issue with l2gmap) 28 - Better management in PCIS code 29 - BDDC with MG framework? 30 31 FETIDP 32 - Move FETIDP code to its own classes 33 34 MATIS related operations contained in BDDC code 35 - Provide general case for subassembling 36 - Preallocation routines in MatISGetMPIAXAIJ 37 38 */ 39 40 /* ---------------------------------------------------------------------------------------------------------------------------------------------- 41 Implementation of BDDC preconditioner based on: 42 C. Dohrmann "An approximate BDDC preconditioner", Numerical Linear Algebra with Applications Volume 14, Issue 2, pages 149-168, March 2007 43 ---------------------------------------------------------------------------------------------------------------------------------------------- */ 44 45 #include "bddc.h" /*I "petscpc.h" I*/ /* includes for fortran wrappers */ 46 #include "bddcprivate.h" 47 #include <petscblaslapack.h> 48 49 /* -------------------------------------------------------------------------- */ 50 #undef __FUNCT__ 51 #define __FUNCT__ "PCSetFromOptions_BDDC" 52 PetscErrorCode PCSetFromOptions_BDDC(PC pc) 53 { 54 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 55 PetscErrorCode ierr; 56 57 PetscFunctionBegin; 58 ierr = PetscOptionsHead("BDDC options");CHKERRQ(ierr); 59 /* Verbose debugging */ 60 ierr = PetscOptionsInt("-pc_bddc_check_level","Verbose output for PCBDDC (intended for debug)","none",pcbddc->dbg_flag,&pcbddc->dbg_flag,NULL);CHKERRQ(ierr); 61 /* Primal space cumstomization */ 62 ierr = PetscOptionsBool("-pc_bddc_use_vertices","Use or not corner dofs in coarse space","none",pcbddc->use_vertices,&pcbddc->use_vertices,NULL);CHKERRQ(ierr); 63 ierr = PetscOptionsBool("-pc_bddc_use_edges","Use or not edge constraints in coarse space","none",pcbddc->use_edges,&pcbddc->use_edges,NULL);CHKERRQ(ierr); 64 ierr = PetscOptionsBool("-pc_bddc_use_faces","Use or not face constraints in coarse space","none",pcbddc->use_faces,&pcbddc->use_faces,NULL);CHKERRQ(ierr); 65 /* Change of basis */ 66 ierr = PetscOptionsBool("-pc_bddc_use_change_of_basis","Use or not change of basis on local edge nodes","none",pcbddc->use_change_of_basis,&pcbddc->use_change_of_basis,NULL);CHKERRQ(ierr); 67 ierr = PetscOptionsBool("-pc_bddc_use_change_on_faces","Use or not change of basis on local face nodes","none",pcbddc->use_change_on_faces,&pcbddc->use_change_on_faces,NULL);CHKERRQ(ierr); 68 if (!pcbddc->use_change_of_basis) { 69 pcbddc->use_change_on_faces = PETSC_FALSE; 70 } 71 /* Switch between M_2 (default) and M_3 preconditioners (as defined by C. Dohrmann in the ref. article) */ 72 ierr = PetscOptionsBool("-pc_bddc_switch_static","Switch on static condensation ops around the interface preconditioner","none",pcbddc->switch_static,&pcbddc->switch_static,NULL);CHKERRQ(ierr); 73 ierr = PetscOptionsInt("-pc_bddc_coarsening_ratio","Set coarsening ratio used in multilevel coarsening","none",pcbddc->coarsening_ratio,&pcbddc->coarsening_ratio,NULL);CHKERRQ(ierr); 74 ierr = PetscOptionsInt("-pc_bddc_levels","Set maximum number of levels for multilevel","none",pcbddc->max_levels,&pcbddc->max_levels,NULL);CHKERRQ(ierr); 75 ierr = PetscOptionsBool("-pc_bddc_use_deluxe_scaling","Use deluxe scaling for BDDC","none",pcbddc->use_deluxe_scaling,&pcbddc->use_deluxe_scaling,NULL);CHKERRQ(ierr); 76 ierr = PetscOptionsTail();CHKERRQ(ierr); 77 PetscFunctionReturn(0); 78 } 79 /* -------------------------------------------------------------------------- */ 80 #undef __FUNCT__ 81 #define __FUNCT__ "PCBDDCSetPrimalVerticesLocalIS_BDDC" 82 static PetscErrorCode PCBDDCSetPrimalVerticesLocalIS_BDDC(PC pc, IS PrimalVertices) 83 { 84 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 85 PetscErrorCode ierr; 86 87 PetscFunctionBegin; 88 ierr = ISDestroy(&pcbddc->user_primal_vertices);CHKERRQ(ierr); 89 ierr = PetscObjectReference((PetscObject)PrimalVertices);CHKERRQ(ierr); 90 pcbddc->user_primal_vertices = PrimalVertices; 91 PetscFunctionReturn(0); 92 } 93 #undef __FUNCT__ 94 #define __FUNCT__ "PCBDDCSetPrimalVerticesLocalIS" 95 /*@ 96 PCBDDCSetPrimalVerticesLocalIS - Set additional user defined primal vertices in PCBDDC 97 98 Not collective 99 100 Input Parameters: 101 + pc - the preconditioning context 102 - PrimalVertices - index set of primal vertices in local numbering 103 104 Level: intermediate 105 106 Notes: 107 108 .seealso: PCBDDC 109 @*/ 110 PetscErrorCode PCBDDCSetPrimalVerticesLocalIS(PC pc, IS PrimalVertices) 111 { 112 PetscErrorCode ierr; 113 114 PetscFunctionBegin; 115 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 116 PetscValidHeaderSpecific(PrimalVertices,IS_CLASSID,2); 117 ierr = PetscTryMethod(pc,"PCBDDCSetPrimalVerticesLocalIS_C",(PC,IS),(pc,PrimalVertices));CHKERRQ(ierr); 118 PetscFunctionReturn(0); 119 } 120 /* -------------------------------------------------------------------------- */ 121 #undef __FUNCT__ 122 #define __FUNCT__ "PCBDDCSetCoarseningRatio_BDDC" 123 static PetscErrorCode PCBDDCSetCoarseningRatio_BDDC(PC pc,PetscInt k) 124 { 125 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 126 127 PetscFunctionBegin; 128 pcbddc->coarsening_ratio = k; 129 PetscFunctionReturn(0); 130 } 131 132 #undef __FUNCT__ 133 #define __FUNCT__ "PCBDDCSetCoarseningRatio" 134 /*@ 135 PCBDDCSetCoarseningRatio - Set coarsening ratio used in multilevel 136 137 Logically collective on PC 138 139 Input Parameters: 140 + pc - the preconditioning context 141 - k - coarsening ratio (H/h at the coarser level) 142 143 Approximatively k subdomains at the finer level will be aggregated into a single subdomain at the coarser level 144 145 Level: intermediate 146 147 Notes: 148 149 .seealso: PCBDDC 150 @*/ 151 PetscErrorCode PCBDDCSetCoarseningRatio(PC pc,PetscInt k) 152 { 153 PetscErrorCode ierr; 154 155 PetscFunctionBegin; 156 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 157 PetscValidLogicalCollectiveInt(pc,k,2); 158 ierr = PetscTryMethod(pc,"PCBDDCSetCoarseningRatio_C",(PC,PetscInt),(pc,k));CHKERRQ(ierr); 159 PetscFunctionReturn(0); 160 } 161 162 /* The following functions (PCBDDCSetUseExactDirichlet PCBDDCSetLevel) are not public */ 163 #undef __FUNCT__ 164 #define __FUNCT__ "PCBDDCSetUseExactDirichlet_BDDC" 165 static PetscErrorCode PCBDDCSetUseExactDirichlet_BDDC(PC pc,PetscBool flg) 166 { 167 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 168 169 PetscFunctionBegin; 170 pcbddc->use_exact_dirichlet_trick = flg; 171 PetscFunctionReturn(0); 172 } 173 174 #undef __FUNCT__ 175 #define __FUNCT__ "PCBDDCSetUseExactDirichlet" 176 PetscErrorCode PCBDDCSetUseExactDirichlet(PC pc,PetscBool flg) 177 { 178 PetscErrorCode ierr; 179 180 PetscFunctionBegin; 181 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 182 PetscValidLogicalCollectiveBool(pc,flg,2); 183 ierr = PetscTryMethod(pc,"PCBDDCSetUseExactDirichlet_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr); 184 PetscFunctionReturn(0); 185 } 186 187 #undef __FUNCT__ 188 #define __FUNCT__ "PCBDDCSetLevel_BDDC" 189 static PetscErrorCode PCBDDCSetLevel_BDDC(PC pc,PetscInt level) 190 { 191 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 192 193 PetscFunctionBegin; 194 pcbddc->current_level = level; 195 PetscFunctionReturn(0); 196 } 197 198 #undef __FUNCT__ 199 #define __FUNCT__ "PCBDDCSetLevel" 200 PetscErrorCode PCBDDCSetLevel(PC pc,PetscInt level) 201 { 202 PetscErrorCode ierr; 203 204 PetscFunctionBegin; 205 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 206 PetscValidLogicalCollectiveInt(pc,level,2); 207 ierr = PetscTryMethod(pc,"PCBDDCSetLevel_C",(PC,PetscInt),(pc,level));CHKERRQ(ierr); 208 PetscFunctionReturn(0); 209 } 210 211 #undef __FUNCT__ 212 #define __FUNCT__ "PCBDDCSetLevels_BDDC" 213 static PetscErrorCode PCBDDCSetLevels_BDDC(PC pc,PetscInt levels) 214 { 215 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 216 217 PetscFunctionBegin; 218 pcbddc->max_levels = levels; 219 PetscFunctionReturn(0); 220 } 221 222 #undef __FUNCT__ 223 #define __FUNCT__ "PCBDDCSetLevels" 224 /*@ 225 PCBDDCSetLevels - Sets the maximum number of levels for multilevel 226 227 Logically collective on PC 228 229 Input Parameters: 230 + pc - the preconditioning context 231 - levels - the maximum number of levels (max 9) 232 233 Default value is 0, i.e. traditional one-level BDDC 234 235 Level: intermediate 236 237 Notes: 238 239 .seealso: PCBDDC 240 @*/ 241 PetscErrorCode PCBDDCSetLevels(PC pc,PetscInt levels) 242 { 243 PetscErrorCode ierr; 244 245 PetscFunctionBegin; 246 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 247 PetscValidLogicalCollectiveInt(pc,levels,2); 248 ierr = PetscTryMethod(pc,"PCBDDCSetLevels_C",(PC,PetscInt),(pc,levels));CHKERRQ(ierr); 249 PetscFunctionReturn(0); 250 } 251 /* -------------------------------------------------------------------------- */ 252 253 #undef __FUNCT__ 254 #define __FUNCT__ "PCBDDCSetNullSpace_BDDC" 255 static PetscErrorCode PCBDDCSetNullSpace_BDDC(PC pc,MatNullSpace NullSpace) 256 { 257 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 258 PetscErrorCode ierr; 259 260 PetscFunctionBegin; 261 ierr = PetscObjectReference((PetscObject)NullSpace);CHKERRQ(ierr); 262 ierr = MatNullSpaceDestroy(&pcbddc->NullSpace);CHKERRQ(ierr); 263 pcbddc->NullSpace = NullSpace; 264 PetscFunctionReturn(0); 265 } 266 267 #undef __FUNCT__ 268 #define __FUNCT__ "PCBDDCSetNullSpace" 269 /*@ 270 PCBDDCSetNullSpace - Set nullspace for BDDC operator 271 272 Logically collective on PC and MatNullSpace 273 274 Input Parameters: 275 + pc - the preconditioning context 276 - NullSpace - Null space of the linear operator to be preconditioned (Pmat) 277 278 Level: intermediate 279 280 Notes: 281 282 .seealso: PCBDDC 283 @*/ 284 PetscErrorCode PCBDDCSetNullSpace(PC pc,MatNullSpace NullSpace) 285 { 286 PetscErrorCode ierr; 287 288 PetscFunctionBegin; 289 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 290 PetscValidHeaderSpecific(NullSpace,MAT_NULLSPACE_CLASSID,2); 291 PetscCheckSameComm(pc,1,NullSpace,2); 292 ierr = PetscTryMethod(pc,"PCBDDCSetNullSpace_C",(PC,MatNullSpace),(pc,NullSpace));CHKERRQ(ierr); 293 PetscFunctionReturn(0); 294 } 295 /* -------------------------------------------------------------------------- */ 296 297 #undef __FUNCT__ 298 #define __FUNCT__ "PCBDDCSetDirichletBoundaries_BDDC" 299 static PetscErrorCode PCBDDCSetDirichletBoundaries_BDDC(PC pc,IS DirichletBoundaries) 300 { 301 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 302 PetscErrorCode ierr; 303 304 PetscFunctionBegin; 305 /* last user setting takes precendence -> destroy any other customization */ 306 ierr = ISDestroy(&pcbddc->DirichletBoundariesLocal);CHKERRQ(ierr); 307 ierr = ISDestroy(&pcbddc->DirichletBoundaries);CHKERRQ(ierr); 308 ierr = PetscObjectReference((PetscObject)DirichletBoundaries);CHKERRQ(ierr); 309 pcbddc->DirichletBoundaries = DirichletBoundaries; 310 pcbddc->recompute_topography = PETSC_TRUE; 311 PetscFunctionReturn(0); 312 } 313 314 #undef __FUNCT__ 315 #define __FUNCT__ "PCBDDCSetDirichletBoundaries" 316 /*@ 317 PCBDDCSetDirichletBoundaries - Set IS defining Dirichlet boundaries for the global problem. 318 319 Collective 320 321 Input Parameters: 322 + pc - the preconditioning context 323 - DirichletBoundaries - parallel IS defining the Dirichlet boundaries 324 325 Level: intermediate 326 327 Notes: Any process can list any global node 328 329 .seealso: PCBDDC 330 @*/ 331 PetscErrorCode PCBDDCSetDirichletBoundaries(PC pc,IS DirichletBoundaries) 332 { 333 PetscErrorCode ierr; 334 335 PetscFunctionBegin; 336 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 337 PetscValidHeaderSpecific(DirichletBoundaries,IS_CLASSID,2); 338 PetscCheckSameComm(pc,1,DirichletBoundaries,2); 339 ierr = PetscTryMethod(pc,"PCBDDCSetDirichletBoundaries_C",(PC,IS),(pc,DirichletBoundaries));CHKERRQ(ierr); 340 PetscFunctionReturn(0); 341 } 342 /* -------------------------------------------------------------------------- */ 343 344 #undef __FUNCT__ 345 #define __FUNCT__ "PCBDDCSetDirichletBoundariesLocal_BDDC" 346 static PetscErrorCode PCBDDCSetDirichletBoundariesLocal_BDDC(PC pc,IS DirichletBoundaries) 347 { 348 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 349 PetscErrorCode ierr; 350 351 PetscFunctionBegin; 352 /* last user setting takes precendence -> destroy any other customization */ 353 ierr = ISDestroy(&pcbddc->DirichletBoundariesLocal);CHKERRQ(ierr); 354 ierr = ISDestroy(&pcbddc->DirichletBoundaries);CHKERRQ(ierr); 355 ierr = PetscObjectReference((PetscObject)DirichletBoundaries);CHKERRQ(ierr); 356 pcbddc->DirichletBoundariesLocal = DirichletBoundaries; 357 PetscFunctionReturn(0); 358 } 359 360 #undef __FUNCT__ 361 #define __FUNCT__ "PCBDDCSetDirichletBoundariesLocal" 362 /*@ 363 PCBDDCSetDirichletBoundariesLocal - Set IS defining Dirichlet boundaries for the global problem in local ordering. 364 365 Collective 366 367 Input Parameters: 368 + pc - the preconditioning context 369 - DirichletBoundaries - parallel IS defining the Dirichlet boundaries (in local ordering) 370 371 Level: intermediate 372 373 Notes: 374 375 .seealso: PCBDDC 376 @*/ 377 PetscErrorCode PCBDDCSetDirichletBoundariesLocal(PC pc,IS DirichletBoundaries) 378 { 379 PetscErrorCode ierr; 380 381 PetscFunctionBegin; 382 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 383 PetscValidHeaderSpecific(DirichletBoundaries,IS_CLASSID,2); 384 PetscCheckSameComm(pc,1,DirichletBoundaries,2); 385 ierr = PetscTryMethod(pc,"PCBDDCSetDirichletBoundariesLocal_C",(PC,IS),(pc,DirichletBoundaries));CHKERRQ(ierr); 386 PetscFunctionReturn(0); 387 } 388 /* -------------------------------------------------------------------------- */ 389 390 #undef __FUNCT__ 391 #define __FUNCT__ "PCBDDCSetNeumannBoundaries_BDDC" 392 static PetscErrorCode PCBDDCSetNeumannBoundaries_BDDC(PC pc,IS NeumannBoundaries) 393 { 394 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 395 PetscErrorCode ierr; 396 397 PetscFunctionBegin; 398 /* last user setting takes precendence -> destroy any other customization */ 399 ierr = ISDestroy(&pcbddc->NeumannBoundariesLocal);CHKERRQ(ierr); 400 ierr = ISDestroy(&pcbddc->NeumannBoundaries);CHKERRQ(ierr); 401 ierr = PetscObjectReference((PetscObject)NeumannBoundaries);CHKERRQ(ierr); 402 pcbddc->NeumannBoundaries = NeumannBoundaries; 403 pcbddc->recompute_topography = PETSC_TRUE; 404 PetscFunctionReturn(0); 405 } 406 407 #undef __FUNCT__ 408 #define __FUNCT__ "PCBDDCSetNeumannBoundaries" 409 /*@ 410 PCBDDCSetNeumannBoundaries - Set IS defining Neumann boundaries for the global problem. 411 412 Collective 413 414 Input Parameters: 415 + pc - the preconditioning context 416 - NeumannBoundaries - parallel IS defining the Neumann boundaries 417 418 Level: intermediate 419 420 Notes: Any process can list any global node 421 422 .seealso: PCBDDC 423 @*/ 424 PetscErrorCode PCBDDCSetNeumannBoundaries(PC pc,IS NeumannBoundaries) 425 { 426 PetscErrorCode ierr; 427 428 PetscFunctionBegin; 429 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 430 PetscValidHeaderSpecific(NeumannBoundaries,IS_CLASSID,2); 431 PetscCheckSameComm(pc,1,NeumannBoundaries,2); 432 ierr = PetscTryMethod(pc,"PCBDDCSetNeumannBoundaries_C",(PC,IS),(pc,NeumannBoundaries));CHKERRQ(ierr); 433 PetscFunctionReturn(0); 434 } 435 /* -------------------------------------------------------------------------- */ 436 437 #undef __FUNCT__ 438 #define __FUNCT__ "PCBDDCSetNeumannBoundariesLocal_BDDC" 439 static PetscErrorCode PCBDDCSetNeumannBoundariesLocal_BDDC(PC pc,IS NeumannBoundaries) 440 { 441 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 442 PetscErrorCode ierr; 443 444 PetscFunctionBegin; 445 /* last user setting takes precendence -> destroy any other customization */ 446 ierr = ISDestroy(&pcbddc->NeumannBoundariesLocal);CHKERRQ(ierr); 447 ierr = ISDestroy(&pcbddc->NeumannBoundaries);CHKERRQ(ierr); 448 ierr = PetscObjectReference((PetscObject)NeumannBoundaries);CHKERRQ(ierr); 449 pcbddc->NeumannBoundariesLocal = NeumannBoundaries; 450 PetscFunctionReturn(0); 451 } 452 453 #undef __FUNCT__ 454 #define __FUNCT__ "PCBDDCSetNeumannBoundariesLocal" 455 /*@ 456 PCBDDCSetNeumannBoundariesLocal - Set IS defining Neumann boundaries for the global problem in local ordering. 457 458 Collective 459 460 Input Parameters: 461 + pc - the preconditioning context 462 - NeumannBoundaries - parallel IS defining the subdomain part of Neumann boundaries (in local ordering) 463 464 Level: intermediate 465 466 Notes: 467 468 .seealso: PCBDDC 469 @*/ 470 PetscErrorCode PCBDDCSetNeumannBoundariesLocal(PC pc,IS NeumannBoundaries) 471 { 472 PetscErrorCode ierr; 473 474 PetscFunctionBegin; 475 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 476 PetscValidHeaderSpecific(NeumannBoundaries,IS_CLASSID,2); 477 PetscCheckSameComm(pc,1,NeumannBoundaries,2); 478 ierr = PetscTryMethod(pc,"PCBDDCSetNeumannBoundariesLocal_C",(PC,IS),(pc,NeumannBoundaries));CHKERRQ(ierr); 479 PetscFunctionReturn(0); 480 } 481 /* -------------------------------------------------------------------------- */ 482 483 #undef __FUNCT__ 484 #define __FUNCT__ "PCBDDCGetDirichletBoundaries_BDDC" 485 static PetscErrorCode PCBDDCGetDirichletBoundaries_BDDC(PC pc,IS *DirichletBoundaries) 486 { 487 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 488 489 PetscFunctionBegin; 490 *DirichletBoundaries = pcbddc->DirichletBoundaries; 491 PetscFunctionReturn(0); 492 } 493 494 #undef __FUNCT__ 495 #define __FUNCT__ "PCBDDCGetDirichletBoundaries" 496 /*@ 497 PCBDDCGetDirichletBoundaries - Get parallel IS for Dirichlet boundaries 498 499 Collective 500 501 Input Parameters: 502 . pc - the preconditioning context 503 504 Output Parameters: 505 . DirichletBoundaries - index set defining the Dirichlet boundaries 506 507 Level: intermediate 508 509 Notes: The IS returned (if any) is the same passed in earlier by the user with PCBDDCSetDirichletBoundaries 510 511 .seealso: PCBDDC 512 @*/ 513 PetscErrorCode PCBDDCGetDirichletBoundaries(PC pc,IS *DirichletBoundaries) 514 { 515 PetscErrorCode ierr; 516 517 PetscFunctionBegin; 518 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 519 ierr = PetscUseMethod(pc,"PCBDDCGetDirichletBoundaries_C",(PC,IS*),(pc,DirichletBoundaries));CHKERRQ(ierr); 520 PetscFunctionReturn(0); 521 } 522 /* -------------------------------------------------------------------------- */ 523 524 #undef __FUNCT__ 525 #define __FUNCT__ "PCBDDCGetDirichletBoundariesLocal_BDDC" 526 static PetscErrorCode PCBDDCGetDirichletBoundariesLocal_BDDC(PC pc,IS *DirichletBoundaries) 527 { 528 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 529 530 PetscFunctionBegin; 531 *DirichletBoundaries = pcbddc->DirichletBoundariesLocal; 532 PetscFunctionReturn(0); 533 } 534 535 #undef __FUNCT__ 536 #define __FUNCT__ "PCBDDCGetDirichletBoundariesLocal" 537 /*@ 538 PCBDDCGetDirichletBoundariesLocal - Get parallel IS for Dirichlet boundaries (in local ordering) 539 540 Collective 541 542 Input Parameters: 543 . pc - the preconditioning context 544 545 Output Parameters: 546 . DirichletBoundaries - index set defining the subdomain part of Dirichlet boundaries 547 548 Level: intermediate 549 550 Notes: 551 552 .seealso: PCBDDC 553 @*/ 554 PetscErrorCode PCBDDCGetDirichletBoundariesLocal(PC pc,IS *DirichletBoundaries) 555 { 556 PetscErrorCode ierr; 557 558 PetscFunctionBegin; 559 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 560 ierr = PetscUseMethod(pc,"PCBDDCGetDirichletBoundariesLocal_C",(PC,IS*),(pc,DirichletBoundaries));CHKERRQ(ierr); 561 PetscFunctionReturn(0); 562 } 563 /* -------------------------------------------------------------------------- */ 564 565 #undef __FUNCT__ 566 #define __FUNCT__ "PCBDDCGetNeumannBoundaries_BDDC" 567 static PetscErrorCode PCBDDCGetNeumannBoundaries_BDDC(PC pc,IS *NeumannBoundaries) 568 { 569 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 570 571 PetscFunctionBegin; 572 *NeumannBoundaries = pcbddc->NeumannBoundaries; 573 PetscFunctionReturn(0); 574 } 575 576 #undef __FUNCT__ 577 #define __FUNCT__ "PCBDDCGetNeumannBoundaries" 578 /*@ 579 PCBDDCGetNeumannBoundaries - Get parallel IS for Neumann boundaries 580 581 Collective 582 583 Input Parameters: 584 . pc - the preconditioning context 585 586 Output Parameters: 587 . NeumannBoundaries - index set defining the Neumann boundaries 588 589 Level: intermediate 590 591 Notes: The IS returned (if any) is the same passed in earlier by the user with PCBDDCSetNeumannBoundaries 592 593 .seealso: PCBDDC 594 @*/ 595 PetscErrorCode PCBDDCGetNeumannBoundaries(PC pc,IS *NeumannBoundaries) 596 { 597 PetscErrorCode ierr; 598 599 PetscFunctionBegin; 600 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 601 ierr = PetscUseMethod(pc,"PCBDDCGetNeumannBoundaries_C",(PC,IS*),(pc,NeumannBoundaries));CHKERRQ(ierr); 602 PetscFunctionReturn(0); 603 } 604 /* -------------------------------------------------------------------------- */ 605 606 #undef __FUNCT__ 607 #define __FUNCT__ "PCBDDCGetNeumannBoundariesLocal_BDDC" 608 static PetscErrorCode PCBDDCGetNeumannBoundariesLocal_BDDC(PC pc,IS *NeumannBoundaries) 609 { 610 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 611 612 PetscFunctionBegin; 613 *NeumannBoundaries = pcbddc->NeumannBoundariesLocal; 614 PetscFunctionReturn(0); 615 } 616 617 #undef __FUNCT__ 618 #define __FUNCT__ "PCBDDCGetNeumannBoundariesLocal" 619 /*@ 620 PCBDDCGetNeumannBoundariesLocal - Get parallel IS for Neumann boundaries (in local ordering) 621 622 Collective 623 624 Input Parameters: 625 . pc - the preconditioning context 626 627 Output Parameters: 628 . NeumannBoundaries - index set defining the subdomain part of Neumann boundaries 629 630 Level: intermediate 631 632 Notes: 633 634 .seealso: PCBDDC 635 @*/ 636 PetscErrorCode PCBDDCGetNeumannBoundariesLocal(PC pc,IS *NeumannBoundaries) 637 { 638 PetscErrorCode ierr; 639 640 PetscFunctionBegin; 641 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 642 ierr = PetscUseMethod(pc,"PCBDDCGetNeumannBoundariesLocal_C",(PC,IS*),(pc,NeumannBoundaries));CHKERRQ(ierr); 643 PetscFunctionReturn(0); 644 } 645 /* -------------------------------------------------------------------------- */ 646 647 #undef __FUNCT__ 648 #define __FUNCT__ "PCBDDCSetLocalAdjacencyGraph_BDDC" 649 static PetscErrorCode PCBDDCSetLocalAdjacencyGraph_BDDC(PC pc, PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode) 650 { 651 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 652 PCBDDCGraph mat_graph = pcbddc->mat_graph; 653 PetscErrorCode ierr; 654 655 PetscFunctionBegin; 656 /* free old CSR */ 657 ierr = PCBDDCGraphResetCSR(mat_graph);CHKERRQ(ierr); 658 /* TODO: PCBDDCGraphSetAdjacency */ 659 /* get CSR into graph structure */ 660 if (copymode == PETSC_COPY_VALUES) { 661 ierr = PetscMalloc((nvtxs+1)*sizeof(PetscInt),&mat_graph->xadj);CHKERRQ(ierr); 662 ierr = PetscMalloc(xadj[nvtxs]*sizeof(PetscInt),&mat_graph->adjncy);CHKERRQ(ierr); 663 ierr = PetscMemcpy(mat_graph->xadj,xadj,(nvtxs+1)*sizeof(PetscInt));CHKERRQ(ierr); 664 ierr = PetscMemcpy(mat_graph->adjncy,adjncy,xadj[nvtxs]*sizeof(PetscInt));CHKERRQ(ierr); 665 } else if (copymode == PETSC_OWN_POINTER) { 666 mat_graph->xadj = (PetscInt*)xadj; 667 mat_graph->adjncy = (PetscInt*)adjncy; 668 } 669 mat_graph->nvtxs_csr = nvtxs; 670 PetscFunctionReturn(0); 671 } 672 673 #undef __FUNCT__ 674 #define __FUNCT__ "PCBDDCSetLocalAdjacencyGraph" 675 /*@ 676 PCBDDCSetLocalAdjacencyGraph - Set adjacency structure (CSR graph) of the local Neumann matrix 677 678 Not collective 679 680 Input Parameters: 681 + pc - the preconditioning context 682 . nvtxs - number of local vertices of the graph (i.e., the local size of your problem) 683 . xadj, adjncy - the CSR graph 684 - copymode - either PETSC_COPY_VALUES or PETSC_OWN_POINTER. 685 686 Level: intermediate 687 688 Notes: 689 690 .seealso: PCBDDC,PetscCopyMode 691 @*/ 692 PetscErrorCode PCBDDCSetLocalAdjacencyGraph(PC pc,PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode) 693 { 694 void (*f)(void) = 0; 695 PetscErrorCode ierr; 696 697 PetscFunctionBegin; 698 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 699 PetscValidIntPointer(xadj,3); 700 PetscValidIntPointer(xadj,4); 701 if (copymode != PETSC_COPY_VALUES && copymode != PETSC_OWN_POINTER) { 702 SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported copy mode %d in %s\n",copymode,__FUNCT__); 703 } 704 ierr = PetscTryMethod(pc,"PCBDDCSetLocalAdjacencyGraph_C",(PC,PetscInt,const PetscInt[],const PetscInt[],PetscCopyMode),(pc,nvtxs,xadj,adjncy,copymode));CHKERRQ(ierr); 705 /* free arrays if PCBDDC is not the PC type */ 706 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",&f);CHKERRQ(ierr); 707 if (!f && copymode == PETSC_OWN_POINTER) { 708 ierr = PetscFree(xadj);CHKERRQ(ierr); 709 ierr = PetscFree(adjncy);CHKERRQ(ierr); 710 } 711 PetscFunctionReturn(0); 712 } 713 /* -------------------------------------------------------------------------- */ 714 715 #undef __FUNCT__ 716 #define __FUNCT__ "PCBDDCSetDofsSplittingLocal_BDDC" 717 static PetscErrorCode PCBDDCSetDofsSplittingLocal_BDDC(PC pc,PetscInt n_is, IS ISForDofs[]) 718 { 719 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 720 PetscInt i; 721 PetscErrorCode ierr; 722 723 PetscFunctionBegin; 724 /* Destroy ISes if they were already set */ 725 for (i=0;i<pcbddc->n_ISForDofsLocal;i++) { 726 ierr = ISDestroy(&pcbddc->ISForDofsLocal[i]);CHKERRQ(ierr); 727 } 728 ierr = PetscFree(pcbddc->ISForDofsLocal);CHKERRQ(ierr); 729 /* last user setting takes precendence -> destroy any other customization */ 730 for (i=0;i<pcbddc->n_ISForDofs;i++) { 731 ierr = ISDestroy(&pcbddc->ISForDofs[i]);CHKERRQ(ierr); 732 } 733 ierr = PetscFree(pcbddc->ISForDofs);CHKERRQ(ierr); 734 pcbddc->n_ISForDofs = 0; 735 /* allocate space then set */ 736 ierr = PetscMalloc(n_is*sizeof(IS),&pcbddc->ISForDofsLocal);CHKERRQ(ierr); 737 for (i=0;i<n_is;i++) { 738 ierr = PetscObjectReference((PetscObject)ISForDofs[i]);CHKERRQ(ierr); 739 pcbddc->ISForDofsLocal[i]=ISForDofs[i]; 740 } 741 pcbddc->n_ISForDofsLocal=n_is; 742 if (n_is) pcbddc->user_provided_isfordofs = PETSC_TRUE; 743 PetscFunctionReturn(0); 744 } 745 746 #undef __FUNCT__ 747 #define __FUNCT__ "PCBDDCSetDofsSplittingLocal" 748 /*@ 749 PCBDDCSetDofsSplittingLocal - Set index sets defining fields of the local subdomain matrix 750 751 Collective 752 753 Input Parameters: 754 + pc - the preconditioning context 755 - n_is - number of index sets defining the fields 756 . ISForDofs - array of IS describing the fields in local ordering 757 758 Level: intermediate 759 760 Notes: n_is should be the same among processes. Not all nodes need to be listed: unlisted nodes will belong to a different field. 761 762 .seealso: PCBDDC 763 @*/ 764 PetscErrorCode PCBDDCSetDofsSplittingLocal(PC pc,PetscInt n_is, IS ISForDofs[]) 765 { 766 PetscInt i; 767 PetscErrorCode ierr; 768 769 PetscFunctionBegin; 770 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 771 PetscValidLogicalCollectiveInt(pc,n_is,2); 772 for (i=0;i<n_is;i++) { 773 PetscCheckSameComm(pc,1,ISForDofs[i],3); 774 PetscValidHeaderSpecific(ISForDofs[i],IS_CLASSID,3); 775 } 776 ierr = PetscTryMethod(pc,"PCBDDCSetDofsSplitting_C",(PC,PetscInt,IS[]),(pc,n_is,ISForDofs));CHKERRQ(ierr); 777 PetscFunctionReturn(0); 778 } 779 /* -------------------------------------------------------------------------- */ 780 781 #undef __FUNCT__ 782 #define __FUNCT__ "PCBDDCSetDofsSplitting_BDDC" 783 static PetscErrorCode PCBDDCSetDofsSplitting_BDDC(PC pc,PetscInt n_is, IS ISForDofs[]) 784 { 785 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 786 PetscInt i; 787 PetscErrorCode ierr; 788 789 PetscFunctionBegin; 790 /* Destroy ISes if they were already set */ 791 for (i=0;i<pcbddc->n_ISForDofs;i++) { 792 ierr = ISDestroy(&pcbddc->ISForDofs[i]);CHKERRQ(ierr); 793 } 794 ierr = PetscFree(pcbddc->ISForDofs);CHKERRQ(ierr); 795 /* last user setting takes precendence -> destroy any other customization */ 796 for (i=0;i<pcbddc->n_ISForDofsLocal;i++) { 797 ierr = ISDestroy(&pcbddc->ISForDofsLocal[i]);CHKERRQ(ierr); 798 } 799 ierr = PetscFree(pcbddc->ISForDofsLocal);CHKERRQ(ierr); 800 pcbddc->n_ISForDofsLocal = 0; 801 /* allocate space then set */ 802 ierr = PetscMalloc(n_is*sizeof(IS),&pcbddc->ISForDofs);CHKERRQ(ierr); 803 for (i=0;i<n_is;i++) { 804 ierr = PetscObjectReference((PetscObject)ISForDofs[i]);CHKERRQ(ierr); 805 pcbddc->ISForDofs[i]=ISForDofs[i]; 806 } 807 pcbddc->n_ISForDofs=n_is; 808 if (n_is) pcbddc->user_provided_isfordofs = PETSC_TRUE; 809 PetscFunctionReturn(0); 810 } 811 812 #undef __FUNCT__ 813 #define __FUNCT__ "PCBDDCSetDofsSplitting" 814 /*@ 815 PCBDDCSetDofsSplitting - Set index sets defining fields of the global matrix 816 817 Collective 818 819 Input Parameters: 820 + pc - the preconditioning context 821 - n_is - number of index sets defining the fields 822 . ISForDofs - array of IS describing the fields in global ordering 823 824 Level: intermediate 825 826 Notes: Any process can list any global node. Not all nodes need to be listed: unlisted nodes will belong to a different field. 827 828 .seealso: PCBDDC 829 @*/ 830 PetscErrorCode PCBDDCSetDofsSplitting(PC pc,PetscInt n_is, IS ISForDofs[]) 831 { 832 PetscInt i; 833 PetscErrorCode ierr; 834 835 PetscFunctionBegin; 836 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 837 PetscValidLogicalCollectiveInt(pc,n_is,2); 838 for (i=0;i<n_is;i++) { 839 PetscCheckSameComm(pc,1,ISForDofs[i],3); 840 PetscValidHeaderSpecific(ISForDofs[i],IS_CLASSID,3); 841 } 842 ierr = PetscTryMethod(pc,"PCBDDCSetDofsSplitting_C",(PC,PetscInt,IS[]),(pc,n_is,ISForDofs));CHKERRQ(ierr); 843 PetscFunctionReturn(0); 844 } 845 /* -------------------------------------------------------------------------- */ 846 #undef __FUNCT__ 847 #define __FUNCT__ "PCPreSolve_BDDC" 848 /* -------------------------------------------------------------------------- */ 849 /* 850 PCPreSolve_BDDC - Changes the right hand side and (if necessary) the initial 851 guess if a transformation of basis approach has been selected. 852 853 Input Parameter: 854 + pc - the preconditioner contex 855 856 Application Interface Routine: PCPreSolve() 857 858 Notes: 859 The interface routine PCPreSolve() is not usually called directly by 860 the user, but instead is called by KSPSolve(). 861 */ 862 static PetscErrorCode PCPreSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x) 863 { 864 PetscErrorCode ierr; 865 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 866 PC_IS *pcis = (PC_IS*)(pc->data); 867 Mat_IS *matis = (Mat_IS*)pc->pmat->data; 868 Mat temp_mat; 869 IS dirIS; 870 Vec used_vec; 871 PetscBool guess_nonzero; 872 873 PetscFunctionBegin; 874 /* if we are working with cg, one dirichlet solve can be avoided during Krylov iterations */ 875 if (ksp) { 876 PetscBool iscg; 877 ierr = PetscObjectTypeCompare((PetscObject)ksp,KSPCG,&iscg);CHKERRQ(ierr); 878 if (!iscg) { 879 ierr = PCBDDCSetUseExactDirichlet(pc,PETSC_FALSE);CHKERRQ(ierr); 880 } 881 } 882 /* Creates parallel work vectors used in presolve */ 883 if (!pcbddc->original_rhs) { 884 ierr = VecDuplicate(pcis->vec1_global,&pcbddc->original_rhs);CHKERRQ(ierr); 885 } 886 if (!pcbddc->temp_solution) { 887 ierr = VecDuplicate(pcis->vec1_global,&pcbddc->temp_solution);CHKERRQ(ierr); 888 } 889 if (x) { 890 ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr); 891 used_vec = x; 892 } else { 893 ierr = PetscObjectReference((PetscObject)pcbddc->temp_solution);CHKERRQ(ierr); 894 used_vec = pcbddc->temp_solution; 895 ierr = VecSet(used_vec,0.0);CHKERRQ(ierr); 896 } 897 /* hack into ksp data structure PCPreSolve comes earlier in src/ksp/ksp/interface/itfunc.c */ 898 if (ksp) { 899 ierr = KSPGetInitialGuessNonzero(ksp,&guess_nonzero);CHKERRQ(ierr); 900 if (!guess_nonzero) { 901 ierr = VecSet(used_vec,0.0);CHKERRQ(ierr); 902 } 903 } 904 905 /* store the original rhs */ 906 ierr = VecCopy(rhs,pcbddc->original_rhs);CHKERRQ(ierr); 907 908 /* Take into account zeroed rows -> change rhs and store solution removed */ 909 /* note that Dirichlet boundaries in global ordering (if any) has already been translated into local ordering in PCBDDCAnalyzeInterface */ 910 ierr = PCBDDCGetDirichletBoundariesLocal(pc,&dirIS);CHKERRQ(ierr); 911 if (rhs && dirIS) { 912 PetscInt dirsize,i,*is_indices; 913 PetscScalar *array_x,*array_diagonal; 914 915 ierr = MatGetDiagonal(pc->pmat,pcis->vec1_global);CHKERRQ(ierr); 916 ierr = VecPointwiseDivide(pcis->vec1_global,rhs,pcis->vec1_global);CHKERRQ(ierr); 917 ierr = VecScatterBegin(matis->ctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 918 ierr = VecScatterEnd(matis->ctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 919 ierr = VecScatterBegin(matis->ctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 920 ierr = VecScatterEnd(matis->ctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 921 ierr = ISGetLocalSize(dirIS,&dirsize);CHKERRQ(ierr); 922 ierr = VecGetArray(pcis->vec1_N,&array_x);CHKERRQ(ierr); 923 ierr = VecGetArray(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr); 924 ierr = ISGetIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr); 925 for (i=0; i<dirsize; i++) array_x[is_indices[i]] = array_diagonal[is_indices[i]]; 926 ierr = ISRestoreIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr); 927 ierr = VecRestoreArray(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr); 928 ierr = VecRestoreArray(pcis->vec1_N,&array_x);CHKERRQ(ierr); 929 ierr = VecScatterBegin(matis->ctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 930 ierr = VecScatterEnd(matis->ctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 931 932 /* remove the computed solution from the rhs */ 933 ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr); 934 ierr = MatMultAdd(pc->pmat,used_vec,rhs,rhs);CHKERRQ(ierr); 935 ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr); 936 } 937 938 /* store partially computed solution and set initial guess */ 939 if (x) { 940 ierr = VecCopy(used_vec,pcbddc->temp_solution);CHKERRQ(ierr); 941 ierr = VecSet(used_vec,0.0);CHKERRQ(ierr); 942 if (pcbddc->use_exact_dirichlet_trick) { 943 ierr = VecScatterBegin(pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 944 ierr = VecScatterEnd (pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 945 ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr); 946 ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 947 ierr = VecScatterEnd (pcis->global_to_D,pcis->vec2_D,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 948 if (ksp) { 949 ierr = KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);CHKERRQ(ierr); 950 } 951 } 952 } 953 954 /* prepare MatMult and rhs for solver */ 955 if (pcbddc->use_change_of_basis) { 956 /* swap pointers for local matrices */ 957 temp_mat = matis->A; 958 matis->A = pcbddc->local_mat; 959 pcbddc->local_mat = temp_mat; 960 if (rhs) { 961 /* Get local rhs and apply transformation of basis */ 962 ierr = VecScatterBegin(pcis->global_to_B,rhs,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 963 ierr = VecScatterEnd (pcis->global_to_B,rhs,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 964 /* from original basis to modified basis */ 965 ierr = MatMultTranspose(pcbddc->ChangeOfBasisMatrix,pcis->vec1_B,pcis->vec2_B);CHKERRQ(ierr); 966 /* put back modified values into the global vec using INSERT_VALUES copy mode */ 967 ierr = VecScatterBegin(pcis->global_to_B,pcis->vec2_B,rhs,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 968 ierr = VecScatterEnd (pcis->global_to_B,pcis->vec2_B,rhs,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 969 } 970 } 971 972 /* remove nullspace if present */ 973 if (ksp && pcbddc->NullSpace) { 974 ierr = MatNullSpaceRemove(pcbddc->NullSpace,used_vec);CHKERRQ(ierr); 975 ierr = MatNullSpaceRemove(pcbddc->NullSpace,rhs);CHKERRQ(ierr); 976 } 977 ierr = VecDestroy(&used_vec);CHKERRQ(ierr); 978 PetscFunctionReturn(0); 979 } 980 /* -------------------------------------------------------------------------- */ 981 #undef __FUNCT__ 982 #define __FUNCT__ "PCPostSolve_BDDC" 983 /* -------------------------------------------------------------------------- */ 984 /* 985 PCPostSolve_BDDC - Changes the computed solution if a transformation of basis 986 approach has been selected. Also, restores rhs to its original state. 987 988 Input Parameter: 989 + pc - the preconditioner contex 990 991 Application Interface Routine: PCPostSolve() 992 993 Notes: 994 The interface routine PCPostSolve() is not usually called directly by 995 the user, but instead is called by KSPSolve(). 996 */ 997 static PetscErrorCode PCPostSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x) 998 { 999 PetscErrorCode ierr; 1000 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 1001 PC_IS *pcis = (PC_IS*)(pc->data); 1002 Mat_IS *matis = (Mat_IS*)pc->pmat->data; 1003 Mat temp_mat; 1004 1005 PetscFunctionBegin; 1006 if (pcbddc->use_change_of_basis) { 1007 /* swap pointers for local matrices */ 1008 temp_mat = matis->A; 1009 matis->A = pcbddc->local_mat; 1010 pcbddc->local_mat = temp_mat; 1011 } 1012 if (pcbddc->use_change_of_basis && x) { 1013 /* Get Local boundary and apply transformation of basis to solution vector */ 1014 ierr = VecScatterBegin(pcis->global_to_B,x,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1015 ierr = VecScatterEnd (pcis->global_to_B,x,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1016 /* from modified basis to original basis */ 1017 ierr = MatMult(pcbddc->ChangeOfBasisMatrix,pcis->vec1_B,pcis->vec2_B);CHKERRQ(ierr); 1018 /* put back modified values into the global vec using INSERT_VALUES copy mode */ 1019 ierr = VecScatterBegin(pcis->global_to_B,pcis->vec2_B,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1020 ierr = VecScatterEnd (pcis->global_to_B,pcis->vec2_B,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1021 } 1022 /* add solution removed in presolve */ 1023 if (x) { 1024 ierr = VecAXPY(x,1.0,pcbddc->temp_solution);CHKERRQ(ierr); 1025 } 1026 /* restore rhs to its original state */ 1027 if (rhs) { 1028 ierr = VecCopy(pcbddc->original_rhs,rhs);CHKERRQ(ierr); 1029 } 1030 PetscFunctionReturn(0); 1031 } 1032 /* -------------------------------------------------------------------------- */ 1033 #undef __FUNCT__ 1034 #define __FUNCT__ "PCSetUp_BDDC" 1035 /* -------------------------------------------------------------------------- */ 1036 /* 1037 PCSetUp_BDDC - Prepares for the use of the BDDC preconditioner 1038 by setting data structures and options. 1039 1040 Input Parameter: 1041 + pc - the preconditioner context 1042 1043 Application Interface Routine: PCSetUp() 1044 1045 Notes: 1046 The interface routine PCSetUp() is not usually called directly by 1047 the user, but instead is called by PCApply() if necessary. 1048 */ 1049 PetscErrorCode PCSetUp_BDDC(PC pc) 1050 { 1051 PetscErrorCode ierr; 1052 PC_BDDC* pcbddc = (PC_BDDC*)pc->data; 1053 MatNullSpace nearnullspace; 1054 MatStructure flag; 1055 PetscBool computeis,computetopography,computesolvers; 1056 PetscBool new_nearnullspace_provided; 1057 1058 PetscFunctionBegin; 1059 /* the following lines of code should be replaced by a better logic between PCIS, PCNN, PCBDDC and other future nonoverlapping preconditioners */ 1060 /* PCIS does not support MatStructure flags different from SAME_PRECONDITIONER */ 1061 /* For BDDC we need to define a local "Neumann" problem different to that defined in PCISSetup 1062 Also, BDDC directly build the Dirichlet problem */ 1063 1064 /* split work */ 1065 if (pc->setupcalled) { 1066 computeis = PETSC_FALSE; 1067 ierr = PCGetOperators(pc,NULL,NULL,&flag);CHKERRQ(ierr); 1068 if (flag == SAME_PRECONDITIONER) { 1069 PetscFunctionReturn(0); 1070 } else if (flag == SAME_NONZERO_PATTERN) { 1071 computetopography = PETSC_FALSE; 1072 computesolvers = PETSC_TRUE; 1073 } else { /* DIFFERENT_NONZERO_PATTERN */ 1074 computetopography = PETSC_TRUE; 1075 computesolvers = PETSC_TRUE; 1076 } 1077 } else { 1078 computeis = PETSC_TRUE; 1079 computetopography = PETSC_TRUE; 1080 computesolvers = PETSC_TRUE; 1081 } 1082 if (pcbddc->recompute_topography) { 1083 computetopography = PETSC_TRUE; 1084 } 1085 1086 /* Get stdout for dbg */ 1087 if (pcbddc->dbg_flag && !pcbddc->dbg_viewer) { 1088 ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)pc),&pcbddc->dbg_viewer);CHKERRQ(ierr); 1089 ierr = PetscViewerASCIISynchronizedAllow(pcbddc->dbg_viewer,PETSC_TRUE);CHKERRQ(ierr); 1090 if (pcbddc->current_level) { 1091 ierr = PetscViewerASCIIAddTab(pcbddc->dbg_viewer,2);CHKERRQ(ierr); 1092 } 1093 } 1094 1095 /* Set up all the "iterative substructuring" common block without computing solvers */ 1096 if (computeis) { 1097 /* HACK INTO PCIS */ 1098 PC_IS* pcis = (PC_IS*)pc->data; 1099 pcis->computesolvers = PETSC_FALSE; 1100 ierr = PCISSetUp(pc);CHKERRQ(ierr); 1101 ierr = ISLocalToGlobalMappingCreateIS(pcis->is_B_local,&pcbddc->BtoNmap);CHKERRQ(ierr); 1102 } 1103 1104 /* Analyze interface */ 1105 if (computetopography) { 1106 ierr = PCBDDCAnalyzeInterface(pc);CHKERRQ(ierr); 1107 } 1108 1109 /* infer if NullSpace object attached to Mat via MatSetNearNullSpace has changed */ 1110 new_nearnullspace_provided = PETSC_FALSE; 1111 ierr = MatGetNearNullSpace(pc->pmat,&nearnullspace);CHKERRQ(ierr); 1112 if (pcbddc->onearnullspace) { /* already used nearnullspace */ 1113 if (!nearnullspace) { /* near null space attached to mat has been destroyed */ 1114 new_nearnullspace_provided = PETSC_TRUE; 1115 } else { 1116 /* determine if the two nullspaces are different (should be lightweight) */ 1117 if (nearnullspace != pcbddc->onearnullspace) { 1118 new_nearnullspace_provided = PETSC_TRUE; 1119 } else { /* maybe the user has changed the content of the nearnullspace so check vectors ObjectStateId */ 1120 PetscInt i; 1121 const Vec *nearnullvecs; 1122 PetscObjectState state; 1123 PetscInt nnsp_size; 1124 ierr = MatNullSpaceGetVecs(nearnullspace,NULL,&nnsp_size,&nearnullvecs);CHKERRQ(ierr); 1125 for (i=0;i<nnsp_size;i++) { 1126 ierr = PetscObjectStateGet((PetscObject)nearnullvecs[i],&state);CHKERRQ(ierr); 1127 if (pcbddc->onearnullvecs_state[i] != state) { 1128 new_nearnullspace_provided = PETSC_TRUE; 1129 break; 1130 } 1131 } 1132 } 1133 } 1134 } else { 1135 if (!nearnullspace) { /* both nearnullspaces are null */ 1136 new_nearnullspace_provided = PETSC_FALSE; 1137 } else { /* nearnullspace attached later */ 1138 new_nearnullspace_provided = PETSC_TRUE; 1139 } 1140 } 1141 1142 /* Setup constraints and related work vectors */ 1143 /* reset primal space flags */ 1144 pcbddc->new_primal_space = PETSC_FALSE; 1145 pcbddc->new_primal_space_local = PETSC_FALSE; 1146 if (computetopography || new_nearnullspace_provided) { 1147 /* It also sets the primal space flags */ 1148 ierr = PCBDDCConstraintsSetUp(pc);CHKERRQ(ierr); 1149 /* Allocate needed local vectors (which depends on quantities defined during ConstraintsSetUp) */ 1150 ierr = PCBDDCSetUpLocalWorkVectors(pc);CHKERRQ(ierr); 1151 } 1152 1153 if (computesolvers || pcbddc->new_primal_space) { 1154 /* reset data */ 1155 ierr = PCBDDCScalingDestroy(pc);CHKERRQ(ierr); 1156 /* Create coarse and local stuffs */ 1157 ierr = PCBDDCSetUpSolvers(pc);CHKERRQ(ierr); 1158 ierr = PCBDDCScalingSetUp(pc);CHKERRQ(ierr); 1159 } 1160 if (pcbddc->dbg_flag && pcbddc->current_level) { 1161 ierr = PetscViewerASCIISubtractTab(pcbddc->dbg_viewer,2);CHKERRQ(ierr); 1162 } 1163 PetscFunctionReturn(0); 1164 } 1165 1166 /* -------------------------------------------------------------------------- */ 1167 /* 1168 PCApply_BDDC - Applies the BDDC preconditioner to a vector. 1169 1170 Input Parameters: 1171 . pc - the preconditioner context 1172 . r - input vector (global) 1173 1174 Output Parameter: 1175 . z - output vector (global) 1176 1177 Application Interface Routine: PCApply() 1178 */ 1179 #undef __FUNCT__ 1180 #define __FUNCT__ "PCApply_BDDC" 1181 PetscErrorCode PCApply_BDDC(PC pc,Vec r,Vec z) 1182 { 1183 PC_IS *pcis = (PC_IS*)(pc->data); 1184 PC_BDDC *pcbddc = (PC_BDDC*)(pc->data); 1185 PetscErrorCode ierr; 1186 const PetscScalar one = 1.0; 1187 const PetscScalar m_one = -1.0; 1188 const PetscScalar zero = 0.0; 1189 1190 /* This code is similar to that provided in nn.c for PCNN 1191 NN interface preconditioner changed to BDDC 1192 Added support for M_3 preconditioner in the reference article (code is active if pcbddc->switch_static = PETSC_TRUE) */ 1193 1194 PetscFunctionBegin; 1195 if (!pcbddc->use_exact_dirichlet_trick) { 1196 /* First Dirichlet solve */ 1197 ierr = VecScatterBegin(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1198 ierr = VecScatterEnd (pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1199 ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr); 1200 /* 1201 Assembling right hand side for BDDC operator 1202 - pcis->vec1_D for the Dirichlet part (if needed, i.e. prec_flag=PETSC_TRUE) 1203 - pcis->vec1_B the interface part of the global vector z 1204 */ 1205 ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr); 1206 ierr = MatMult(pcis->A_BI,pcis->vec2_D,pcis->vec1_B);CHKERRQ(ierr); 1207 if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec2_D,pcis->vec1_D,pcis->vec1_D);CHKERRQ(ierr); } 1208 ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr); 1209 ierr = VecCopy(r,z);CHKERRQ(ierr); 1210 ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1211 ierr = VecScatterEnd (pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1212 ierr = PCBDDCScalingRestriction(pc,z,pcis->vec1_B);CHKERRQ(ierr); 1213 } else { 1214 ierr = VecSet(pcis->vec1_D,zero);CHKERRQ(ierr); 1215 ierr = VecSet(pcis->vec2_D,zero);CHKERRQ(ierr); 1216 ierr = PCBDDCScalingRestriction(pc,r,pcis->vec1_B);CHKERRQ(ierr); 1217 } 1218 1219 /* Apply interface preconditioner 1220 input/output vecs: pcis->vec1_B and pcis->vec1_D */ 1221 ierr = PCBDDCApplyInterfacePreconditioner(pc);CHKERRQ(ierr); 1222 1223 /* Apply transpose of partition of unity operator */ 1224 ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,z);CHKERRQ(ierr); 1225 1226 /* Second Dirichlet solve and assembling of output */ 1227 ierr = VecScatterBegin(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1228 ierr = VecScatterEnd (pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1229 ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec3_D);CHKERRQ(ierr); 1230 if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec1_D,pcis->vec3_D,pcis->vec3_D);CHKERRQ(ierr); } 1231 ierr = KSPSolve(pcbddc->ksp_D,pcis->vec3_D,pcis->vec4_D);CHKERRQ(ierr); 1232 ierr = VecScale(pcis->vec4_D,m_one);CHKERRQ(ierr); 1233 if (pcbddc->switch_static) { ierr = VecAXPY (pcis->vec4_D,one,pcis->vec1_D);CHKERRQ(ierr); } 1234 ierr = VecAXPY (pcis->vec2_D,one,pcis->vec4_D);CHKERRQ(ierr); 1235 ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1236 ierr = VecScatterEnd (pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1237 PetscFunctionReturn(0); 1238 } 1239 /* -------------------------------------------------------------------------- */ 1240 1241 #undef __FUNCT__ 1242 #define __FUNCT__ "PCDestroy_BDDC" 1243 PetscErrorCode PCDestroy_BDDC(PC pc) 1244 { 1245 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 1246 PetscErrorCode ierr; 1247 1248 PetscFunctionBegin; 1249 /* free data created by PCIS */ 1250 ierr = PCISDestroy(pc);CHKERRQ(ierr); 1251 /* free BDDC custom data */ 1252 ierr = PCBDDCResetCustomization(pc);CHKERRQ(ierr); 1253 /* destroy objects related to topography */ 1254 ierr = PCBDDCResetTopography(pc);CHKERRQ(ierr); 1255 /* free allocated graph structure */ 1256 ierr = PetscFree(pcbddc->mat_graph);CHKERRQ(ierr); 1257 /* free data for scaling operator */ 1258 ierr = PCBDDCScalingDestroy(pc);CHKERRQ(ierr); 1259 /* free solvers stuff */ 1260 ierr = PCBDDCResetSolvers(pc);CHKERRQ(ierr); 1261 /* free global vectors needed in presolve */ 1262 ierr = VecDestroy(&pcbddc->temp_solution);CHKERRQ(ierr); 1263 ierr = VecDestroy(&pcbddc->original_rhs);CHKERRQ(ierr); 1264 ierr = ISLocalToGlobalMappingDestroy(&pcbddc->BtoNmap);CHKERRQ(ierr); 1265 /* remove functions */ 1266 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",NULL);CHKERRQ(ierr); 1267 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",NULL);CHKERRQ(ierr); 1268 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",NULL);CHKERRQ(ierr); 1269 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",NULL);CHKERRQ(ierr); 1270 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",NULL);CHKERRQ(ierr); 1271 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",NULL);CHKERRQ(ierr); 1272 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",NULL);CHKERRQ(ierr); 1273 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundariesLocal_C",NULL);CHKERRQ(ierr); 1274 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",NULL);CHKERRQ(ierr); 1275 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundariesLocal_C",NULL);CHKERRQ(ierr); 1276 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",NULL);CHKERRQ(ierr); 1277 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundariesLocal_C",NULL);CHKERRQ(ierr); 1278 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",NULL);CHKERRQ(ierr); 1279 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundariesLocal_C",NULL);CHKERRQ(ierr); 1280 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",NULL);CHKERRQ(ierr); 1281 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplittingLocal_C",NULL);CHKERRQ(ierr); 1282 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",NULL);CHKERRQ(ierr); 1283 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",NULL);CHKERRQ(ierr); 1284 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",NULL);CHKERRQ(ierr); 1285 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",NULL);CHKERRQ(ierr); 1286 /* Free the private data structure */ 1287 ierr = PetscFree(pc->data);CHKERRQ(ierr); 1288 PetscFunctionReturn(0); 1289 } 1290 /* -------------------------------------------------------------------------- */ 1291 1292 #undef __FUNCT__ 1293 #define __FUNCT__ "PCBDDCMatFETIDPGetRHS_BDDC" 1294 static PetscErrorCode PCBDDCMatFETIDPGetRHS_BDDC(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs) 1295 { 1296 FETIDPMat_ctx mat_ctx; 1297 PC_IS* pcis; 1298 PC_BDDC* pcbddc; 1299 PetscErrorCode ierr; 1300 1301 PetscFunctionBegin; 1302 ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr); 1303 pcis = (PC_IS*)mat_ctx->pc->data; 1304 pcbddc = (PC_BDDC*)mat_ctx->pc->data; 1305 1306 /* change of basis for physical rhs if needed 1307 It also changes the rhs in case of dirichlet boundaries */ 1308 ierr = PCPreSolve_BDDC(mat_ctx->pc,NULL,standard_rhs,NULL);CHKERRQ(ierr); 1309 /* store vectors for computation of fetidp final solution */ 1310 ierr = VecScatterBegin(pcis->global_to_D,standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1311 ierr = VecScatterEnd(pcis->global_to_D,standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1312 /* scale rhs since it should be unassembled */ 1313 /* TODO use counter scaling? (also below) */ 1314 ierr = VecScatterBegin(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1315 ierr = VecScatterEnd(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1316 /* Apply partition of unity */ 1317 ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr); 1318 /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */ 1319 if (!pcbddc->switch_static) { 1320 /* compute partially subassembled Schur complement right-hand side */ 1321 ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr); 1322 ierr = MatMult(pcis->A_BI,pcis->vec1_D,pcis->vec1_B);CHKERRQ(ierr); 1323 ierr = VecAXPY(mat_ctx->temp_solution_B,-1.0,pcis->vec1_B);CHKERRQ(ierr); 1324 ierr = VecSet(standard_rhs,0.0);CHKERRQ(ierr); 1325 ierr = VecScatterBegin(pcis->global_to_B,mat_ctx->temp_solution_B,standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1326 ierr = VecScatterEnd(pcis->global_to_B,mat_ctx->temp_solution_B,standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1327 /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */ 1328 ierr = VecScatterBegin(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1329 ierr = VecScatterEnd(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1330 ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr); 1331 } 1332 /* BDDC rhs */ 1333 ierr = VecCopy(mat_ctx->temp_solution_B,pcis->vec1_B);CHKERRQ(ierr); 1334 if (pcbddc->switch_static) { 1335 ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr); 1336 } 1337 /* apply BDDC */ 1338 ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc);CHKERRQ(ierr); 1339 /* Application of B_delta and assembling of rhs for fetidp fluxes */ 1340 ierr = VecSet(fetidp_flux_rhs,0.0);CHKERRQ(ierr); 1341 ierr = MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local);CHKERRQ(ierr); 1342 ierr = VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1343 ierr = VecScatterEnd (mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1344 /* restore original rhs */ 1345 ierr = VecCopy(pcbddc->original_rhs,standard_rhs);CHKERRQ(ierr); 1346 PetscFunctionReturn(0); 1347 } 1348 1349 #undef __FUNCT__ 1350 #define __FUNCT__ "PCBDDCMatFETIDPGetRHS" 1351 /*@ 1352 PCBDDCMatFETIDPGetRHS - Compute the right-hand side for FETIDP linear system 1353 1354 Collective 1355 1356 Input Parameters: 1357 + fetidp_mat - the FETIDP matrix object obtained by calling PCBDDCCreateFETIDPOperators 1358 . standard_rhs - the right-hand side for your linear system 1359 1360 Output Parameters: 1361 - fetidp_flux_rhs - the right-hand side for the FETIDP linear system 1362 1363 Level: developer 1364 1365 Notes: 1366 1367 .seealso: PCBDDC,PCBDDCCreateFETIDPOperators 1368 @*/ 1369 PetscErrorCode PCBDDCMatFETIDPGetRHS(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs) 1370 { 1371 FETIDPMat_ctx mat_ctx; 1372 PetscErrorCode ierr; 1373 1374 PetscFunctionBegin; 1375 ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr); 1376 ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetRHS_C",(Mat,Vec,Vec),(fetidp_mat,standard_rhs,fetidp_flux_rhs));CHKERRQ(ierr); 1377 PetscFunctionReturn(0); 1378 } 1379 /* -------------------------------------------------------------------------- */ 1380 1381 #undef __FUNCT__ 1382 #define __FUNCT__ "PCBDDCMatFETIDPGetSolution_BDDC" 1383 static PetscErrorCode PCBDDCMatFETIDPGetSolution_BDDC(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol) 1384 { 1385 FETIDPMat_ctx mat_ctx; 1386 PC_IS* pcis; 1387 PC_BDDC* pcbddc; 1388 PetscErrorCode ierr; 1389 1390 PetscFunctionBegin; 1391 ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr); 1392 pcis = (PC_IS*)mat_ctx->pc->data; 1393 pcbddc = (PC_BDDC*)mat_ctx->pc->data; 1394 1395 /* apply B_delta^T */ 1396 ierr = VecScatterBegin(mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1397 ierr = VecScatterEnd (mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1398 ierr = MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr); 1399 /* compute rhs for BDDC application */ 1400 ierr = VecAYPX(pcis->vec1_B,-1.0,mat_ctx->temp_solution_B);CHKERRQ(ierr); 1401 if (pcbddc->switch_static) { 1402 ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr); 1403 } 1404 /* apply BDDC */ 1405 ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc);CHKERRQ(ierr); 1406 /* put values into standard global vector */ 1407 ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1408 ierr = VecScatterEnd (pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1409 if (!pcbddc->switch_static) { 1410 /* compute values into the interior if solved for the partially subassembled Schur complement */ 1411 ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec1_D);CHKERRQ(ierr); 1412 ierr = VecAXPY(mat_ctx->temp_solution_D,-1.0,pcis->vec1_D);CHKERRQ(ierr); 1413 ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr); 1414 } 1415 ierr = VecScatterBegin(pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1416 ierr = VecScatterEnd (pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1417 /* final change of basis if needed 1418 Is also sums the dirichlet part removed during RHS assembling */ 1419 ierr = PCPostSolve_BDDC(mat_ctx->pc,NULL,NULL,standard_sol);CHKERRQ(ierr); 1420 PetscFunctionReturn(0); 1421 } 1422 1423 #undef __FUNCT__ 1424 #define __FUNCT__ "PCBDDCMatFETIDPGetSolution" 1425 /*@ 1426 PCBDDCMatFETIDPGetSolution - Compute the physical solution from the solution of the FETIDP linear system 1427 1428 Collective 1429 1430 Input Parameters: 1431 + fetidp_mat - the FETIDP matrix obtained by calling PCBDDCCreateFETIDPOperators 1432 . fetidp_flux_sol - the solution of the FETIDP linear system 1433 1434 Output Parameters: 1435 - standard_sol - the solution defined on the physical domain 1436 1437 Level: developer 1438 1439 Notes: 1440 1441 .seealso: PCBDDC,PCBDDCCreateFETIDPOperators 1442 @*/ 1443 PetscErrorCode PCBDDCMatFETIDPGetSolution(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol) 1444 { 1445 FETIDPMat_ctx mat_ctx; 1446 PetscErrorCode ierr; 1447 1448 PetscFunctionBegin; 1449 ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr); 1450 ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetSolution_C",(Mat,Vec,Vec),(fetidp_mat,fetidp_flux_sol,standard_sol));CHKERRQ(ierr); 1451 PetscFunctionReturn(0); 1452 } 1453 /* -------------------------------------------------------------------------- */ 1454 1455 extern PetscErrorCode FETIDPMatMult(Mat,Vec,Vec); 1456 extern PetscErrorCode PCBDDCDestroyFETIDPMat(Mat); 1457 extern PetscErrorCode FETIDPPCApply(PC,Vec,Vec); 1458 extern PetscErrorCode PCBDDCDestroyFETIDPPC(PC); 1459 1460 #undef __FUNCT__ 1461 #define __FUNCT__ "PCBDDCCreateFETIDPOperators_BDDC" 1462 static PetscErrorCode PCBDDCCreateFETIDPOperators_BDDC(PC pc, Mat *fetidp_mat, PC *fetidp_pc) 1463 { 1464 1465 FETIDPMat_ctx fetidpmat_ctx; 1466 Mat newmat; 1467 FETIDPPC_ctx fetidppc_ctx; 1468 PC newpc; 1469 MPI_Comm comm; 1470 PetscErrorCode ierr; 1471 1472 PetscFunctionBegin; 1473 ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 1474 /* FETIDP linear matrix */ 1475 ierr = PCBDDCCreateFETIDPMatContext(pc,&fetidpmat_ctx);CHKERRQ(ierr); 1476 ierr = PCBDDCSetupFETIDPMatContext(fetidpmat_ctx);CHKERRQ(ierr); 1477 ierr = MatCreateShell(comm,PETSC_DECIDE,PETSC_DECIDE,fetidpmat_ctx->n_lambda,fetidpmat_ctx->n_lambda,fetidpmat_ctx,&newmat);CHKERRQ(ierr); 1478 ierr = MatShellSetOperation(newmat,MATOP_MULT,(void (*)(void))FETIDPMatMult);CHKERRQ(ierr); 1479 ierr = MatShellSetOperation(newmat,MATOP_DESTROY,(void (*)(void))PCBDDCDestroyFETIDPMat);CHKERRQ(ierr); 1480 ierr = MatSetUp(newmat);CHKERRQ(ierr); 1481 /* FETIDP preconditioner */ 1482 ierr = PCBDDCCreateFETIDPPCContext(pc,&fetidppc_ctx);CHKERRQ(ierr); 1483 ierr = PCBDDCSetupFETIDPPCContext(newmat,fetidppc_ctx);CHKERRQ(ierr); 1484 ierr = PCCreate(comm,&newpc);CHKERRQ(ierr); 1485 ierr = PCSetType(newpc,PCSHELL);CHKERRQ(ierr); 1486 ierr = PCShellSetContext(newpc,fetidppc_ctx);CHKERRQ(ierr); 1487 ierr = PCShellSetApply(newpc,FETIDPPCApply);CHKERRQ(ierr); 1488 ierr = PCShellSetDestroy(newpc,PCBDDCDestroyFETIDPPC);CHKERRQ(ierr); 1489 ierr = PCSetOperators(newpc,newmat,newmat,SAME_PRECONDITIONER);CHKERRQ(ierr); 1490 ierr = PCSetUp(newpc);CHKERRQ(ierr); 1491 /* return pointers for objects created */ 1492 *fetidp_mat=newmat; 1493 *fetidp_pc=newpc; 1494 PetscFunctionReturn(0); 1495 } 1496 1497 #undef __FUNCT__ 1498 #define __FUNCT__ "PCBDDCCreateFETIDPOperators" 1499 /*@ 1500 PCBDDCCreateFETIDPOperators - Create operators for FETIDP 1501 1502 Collective 1503 1504 Input Parameters: 1505 + pc - the BDDC preconditioning context already setup 1506 1507 Output Parameters: 1508 . fetidp_mat - shell FETIDP matrix object 1509 . fetidp_pc - shell Dirichlet preconditioner for FETIDP matrix 1510 1511 Options Database Keys: 1512 - -fetidp_fullyredundant: use or not a fully redundant set of Lagrange multipliers 1513 1514 Level: developer 1515 1516 Notes: 1517 Currently the only operation provided for FETIDP matrix is MatMult 1518 1519 .seealso: PCBDDC 1520 @*/ 1521 PetscErrorCode PCBDDCCreateFETIDPOperators(PC pc, Mat *fetidp_mat, PC *fetidp_pc) 1522 { 1523 PetscErrorCode ierr; 1524 1525 PetscFunctionBegin; 1526 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1527 if (pc->setupcalled) { 1528 ierr = PetscUseMethod(pc,"PCBDDCCreateFETIDPOperators_C",(PC,Mat*,PC*),(pc,fetidp_mat,fetidp_pc));CHKERRQ(ierr); 1529 } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"You must call PCSetup_BDDC() first \n"); 1530 PetscFunctionReturn(0); 1531 } 1532 /* -------------------------------------------------------------------------- */ 1533 /*MC 1534 PCBDDC - Balancing Domain Decomposition by Constraints. 1535 1536 An implementation of the BDDC preconditioner based on 1537 1538 .vb 1539 [1] C. R. Dohrmann. "An approximate BDDC preconditioner", Numerical Linear Algebra with Applications Volume 14, Issue 2, pages 149-168, March 2007 1540 [2] A. Klawonn and O. B. Widlund. "Dual-Primal FETI Methods for Linear Elasticity", http://cs.nyu.edu/csweb/Research/TechReports/TR2004-855/TR2004-855.pdf 1541 [3] J. Mandel, B. Sousedik, C. R. Dohrmann. "Multispace and Multilevel BDDC", http://arxiv.org/abs/0712.3977 1542 .ve 1543 1544 The matrix to be preconditioned (Pmat) must be of type MATIS. 1545 1546 Currently works with MATIS matrices with local Neumann matrices of type MATSEQAIJ, MATSEQBAIJ or MATSEQSBAIJ, either with real or complex numbers. 1547 1548 It also works with unsymmetric and indefinite problems. 1549 1550 Unlike 'conventional' interface preconditioners, PCBDDC iterates over all degrees of freedom, not just those on the interface. This allows the use of approximate solvers on the subdomains. 1551 1552 Approximate local solvers are automatically adapted for singular linear problems (see [1]) if the user has provided the nullspace using PCBDDCSetNullSpace 1553 1554 Boundary nodes are split in vertices, edges and faces using information from the local to global mapping of dofs and the local connectivity graph of nodes. The latter can be customized by using PCBDDCSetLocalAdjacencyGraph 1555 1556 Constraints can be customized by attaching a MatNullSpace object to the MATIS matrix via MatSetNearNullSpace. 1557 1558 Change of basis is performed similarly to [2] when requested. When more the one constraint is present on a single connected component (i.e. an edge or a face), a robust method based on local QR factorizations is used. 1559 1560 The PETSc implementation also supports multilevel BDDC [3]. Coarse grids are partitioned using MatPartitioning object. 1561 1562 Options Database Keys: 1563 1564 . -pc_bddc_use_vertices <1> - use or not vertices in primal space 1565 . -pc_bddc_use_edges <1> - use or not edges in primal space 1566 . -pc_bddc_use_faces <0> - use or not faces in primal space 1567 . -pc_bddc_use_change_of_basis <0> - use change of basis approach (on edges only) 1568 . -pc_bddc_use_change_on_faces <0> - use change of basis approach on faces if change of basis has been requested 1569 . -pc_bddc_switch_static <0> - switches from M_2 to M_3 operator (see reference article [1]) 1570 . -pc_bddc_levels <0> - maximum number of levels for multilevel 1571 . -pc_bddc_coarsening_ratio - H/h ratio at the coarser level 1572 - -pc_bddc_check_level <0> - set verbosity level of debugging output 1573 1574 Options for Dirichlet, Neumann or coarse solver can be set with 1575 .vb 1576 -pc_bddc_dirichlet_ 1577 -pc_bddc_neumann_ 1578 -pc_bddc_coarse_ 1579 .ve 1580 e.g -pc_bddc_dirichlet_ksp_type richardson -pc_bddc_dirichlet_pc_type gamg 1581 1582 When using a multilevel approach, solvers' options at the N-th level can be specified as 1583 .vb 1584 -pc_bddc_dirichlet_N_ 1585 -pc_bddc_neumann_N_ 1586 -pc_bddc_coarse_N_ 1587 .ve 1588 Note that level number ranges from the finest 0 to the coarsest N 1589 1590 Level: intermediate 1591 1592 Developer notes: 1593 Currently does not work with KSPBCGS and other KSPs requiring the specialization of PCApplyTranspose 1594 1595 New deluxe scaling operator will be available soon. 1596 1597 Contributed by Stefano Zampini 1598 1599 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, MATIS 1600 M*/ 1601 1602 #undef __FUNCT__ 1603 #define __FUNCT__ "PCCreate_BDDC" 1604 PETSC_EXTERN PetscErrorCode PCCreate_BDDC(PC pc) 1605 { 1606 PetscErrorCode ierr; 1607 PC_BDDC *pcbddc; 1608 1609 PetscFunctionBegin; 1610 /* Creates the private data structure for this preconditioner and attach it to the PC object. */ 1611 ierr = PetscNewLog(pc,PC_BDDC,&pcbddc);CHKERRQ(ierr); 1612 pc->data = (void*)pcbddc; 1613 1614 /* create PCIS data structure */ 1615 ierr = PCISCreate(pc);CHKERRQ(ierr); 1616 1617 /* BDDC customization */ 1618 pcbddc->use_vertices = PETSC_TRUE; 1619 pcbddc->use_edges = PETSC_TRUE; 1620 pcbddc->use_faces = PETSC_FALSE; 1621 pcbddc->use_change_of_basis = PETSC_FALSE; 1622 pcbddc->use_change_on_faces = PETSC_FALSE; 1623 pcbddc->switch_static = PETSC_FALSE; 1624 pcbddc->use_nnsp_true = PETSC_FALSE; /* not yet exposed */ 1625 pcbddc->dbg_flag = 0; 1626 1627 pcbddc->BtoNmap = 0; 1628 pcbddc->local_primal_size = 0; 1629 pcbddc->n_vertices = 0; 1630 pcbddc->n_actual_vertices = 0; 1631 pcbddc->n_constraints = 0; 1632 pcbddc->primal_indices_local_idxs = 0; 1633 pcbddc->recompute_topography = PETSC_FALSE; 1634 pcbddc->coarse_size = 0; 1635 pcbddc->new_primal_space = PETSC_FALSE; 1636 pcbddc->new_primal_space_local = PETSC_FALSE; 1637 pcbddc->global_primal_indices = 0; 1638 pcbddc->onearnullspace = 0; 1639 pcbddc->onearnullvecs_state = 0; 1640 pcbddc->user_primal_vertices = 0; 1641 pcbddc->NullSpace = 0; 1642 pcbddc->temp_solution = 0; 1643 pcbddc->original_rhs = 0; 1644 pcbddc->local_mat = 0; 1645 pcbddc->ChangeOfBasisMatrix = 0; 1646 pcbddc->coarse_vec = 0; 1647 pcbddc->coarse_rhs = 0; 1648 pcbddc->coarse_ksp = 0; 1649 pcbddc->coarse_phi_B = 0; 1650 pcbddc->coarse_phi_D = 0; 1651 pcbddc->coarse_psi_B = 0; 1652 pcbddc->coarse_psi_D = 0; 1653 pcbddc->vec1_P = 0; 1654 pcbddc->vec1_R = 0; 1655 pcbddc->vec2_R = 0; 1656 pcbddc->local_auxmat1 = 0; 1657 pcbddc->local_auxmat2 = 0; 1658 pcbddc->R_to_B = 0; 1659 pcbddc->R_to_D = 0; 1660 pcbddc->ksp_D = 0; 1661 pcbddc->ksp_R = 0; 1662 pcbddc->NeumannBoundaries = 0; 1663 pcbddc->NeumannBoundariesLocal = 0; 1664 pcbddc->DirichletBoundaries = 0; 1665 pcbddc->DirichletBoundariesLocal = 0; 1666 pcbddc->user_provided_isfordofs = PETSC_FALSE; 1667 pcbddc->n_ISForDofs = 0; 1668 pcbddc->n_ISForDofsLocal = 0; 1669 pcbddc->ISForDofs = 0; 1670 pcbddc->ISForDofsLocal = 0; 1671 pcbddc->ConstraintMatrix = 0; 1672 pcbddc->use_exact_dirichlet_trick = PETSC_TRUE; 1673 pcbddc->coarse_loc_to_glob = 0; 1674 pcbddc->coarsening_ratio = 8; 1675 pcbddc->current_level = 0; 1676 pcbddc->max_levels = 0; 1677 1678 /* create local graph structure */ 1679 ierr = PCBDDCGraphCreate(&pcbddc->mat_graph);CHKERRQ(ierr); 1680 1681 /* scaling */ 1682 pcbddc->use_deluxe_scaling = PETSC_FALSE; 1683 pcbddc->work_scaling = 0; 1684 1685 /* function pointers */ 1686 pc->ops->apply = PCApply_BDDC; 1687 pc->ops->applytranspose = 0; 1688 pc->ops->setup = PCSetUp_BDDC; 1689 pc->ops->destroy = PCDestroy_BDDC; 1690 pc->ops->setfromoptions = PCSetFromOptions_BDDC; 1691 pc->ops->view = 0; 1692 pc->ops->applyrichardson = 0; 1693 pc->ops->applysymmetricleft = 0; 1694 pc->ops->applysymmetricright = 0; 1695 pc->ops->presolve = PCPreSolve_BDDC; 1696 pc->ops->postsolve = PCPostSolve_BDDC; 1697 1698 /* composing function */ 1699 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",PCBDDCSetPrimalVerticesLocalIS_BDDC);CHKERRQ(ierr); 1700 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",PCBDDCSetCoarseningRatio_BDDC);CHKERRQ(ierr); 1701 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",PCBDDCSetLevel_BDDC);CHKERRQ(ierr); 1702 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",PCBDDCSetUseExactDirichlet_BDDC);CHKERRQ(ierr); 1703 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",PCBDDCSetLevels_BDDC);CHKERRQ(ierr); 1704 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",PCBDDCSetNullSpace_BDDC);CHKERRQ(ierr); 1705 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",PCBDDCSetDirichletBoundaries_BDDC);CHKERRQ(ierr); 1706 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundariesLocal_C",PCBDDCSetDirichletBoundariesLocal_BDDC);CHKERRQ(ierr); 1707 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",PCBDDCSetNeumannBoundaries_BDDC);CHKERRQ(ierr); 1708 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundariesLocal_C",PCBDDCSetNeumannBoundariesLocal_BDDC);CHKERRQ(ierr); 1709 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",PCBDDCGetDirichletBoundaries_BDDC);CHKERRQ(ierr); 1710 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundariesLocal_C",PCBDDCGetDirichletBoundariesLocal_BDDC);CHKERRQ(ierr); 1711 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",PCBDDCGetNeumannBoundaries_BDDC);CHKERRQ(ierr); 1712 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundariesLocal_C",PCBDDCGetNeumannBoundariesLocal_BDDC);CHKERRQ(ierr); 1713 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",PCBDDCSetDofsSplitting_BDDC);CHKERRQ(ierr); 1714 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplittingLocal_C",PCBDDCSetDofsSplittingLocal_BDDC);CHKERRQ(ierr); 1715 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",PCBDDCSetLocalAdjacencyGraph_BDDC);CHKERRQ(ierr); 1716 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",PCBDDCCreateFETIDPOperators_BDDC);CHKERRQ(ierr); 1717 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",PCBDDCMatFETIDPGetRHS_BDDC);CHKERRQ(ierr); 1718 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",PCBDDCMatFETIDPGetSolution_BDDC);CHKERRQ(ierr); 1719 PetscFunctionReturn(0); 1720 } 1721 1722