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