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