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