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