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