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