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