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