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