1 /* TODOLIST 2 Don't use qr when number of primal dof per cc is 1 3 Better management for BAIJ local mats. How to deal with SBAIJ? 4 Provide PCApplyTranpose 5 make runexe59 6 Man pages 7 Move FETIDP code 8 Provide general case for subassembling 9 Preallocation routines in MatConvert_IS_AIJ 10 Allow different customizations between solves 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 user defined primal vertices in PCBDDC. 76 77 Not collective 78 79 Input Parameters: 80 + pc - the preconditioning context 81 - PrimalVertices - index sets 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 coarsening 115 116 Logically collective on PC 117 118 Input Parameters: 119 + pc - the preconditioning context 120 - k - coarsening ratio 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 within the multilevel approach. 205 206 Logically collective on PC 207 208 Input Parameters: 209 + pc - the preconditioning context 210 - levels - the maximum number of levels 211 212 Default value is 0, i.e. coarse problem will be solved exactly 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 of global operator of BDDC preconditioned mat. 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. 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 index set defining subdomain part (in local ordering) 294 of Dirichlet boundaries for the global problem. 295 296 Not collective 297 298 Input Parameters: 299 + pc - the preconditioning context 300 - DirichletBoundaries - sequential index set defining the subdomain part of Dirichlet boundaries (can be NULL) 301 302 Level: intermediate 303 304 Notes: 305 306 .seealso: PCBDDC 307 @*/ 308 PetscErrorCode PCBDDCSetDirichletBoundaries(PC pc,IS DirichletBoundaries) 309 { 310 PetscErrorCode ierr; 311 312 PetscFunctionBegin; 313 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 314 PetscValidHeaderSpecific(DirichletBoundaries,IS_CLASSID,2); 315 ierr = PetscTryMethod(pc,"PCBDDCSetDirichletBoundaries_C",(PC,IS),(pc,DirichletBoundaries));CHKERRQ(ierr); 316 PetscFunctionReturn(0); 317 } 318 /* -------------------------------------------------------------------------- */ 319 320 #undef __FUNCT__ 321 #define __FUNCT__ "PCBDDCSetNeumannBoundaries_BDDC" 322 static PetscErrorCode PCBDDCSetNeumannBoundaries_BDDC(PC pc,IS NeumannBoundaries) 323 { 324 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 325 PetscErrorCode ierr; 326 327 PetscFunctionBegin; 328 ierr = ISDestroy(&pcbddc->NeumannBoundaries);CHKERRQ(ierr); 329 ierr = PetscObjectReference((PetscObject)NeumannBoundaries);CHKERRQ(ierr); 330 pcbddc->NeumannBoundaries=NeumannBoundaries; 331 PetscFunctionReturn(0); 332 } 333 334 #undef __FUNCT__ 335 #define __FUNCT__ "PCBDDCSetNeumannBoundaries" 336 /*@ 337 PCBDDCSetNeumannBoundaries - Set index set defining subdomain part (in local ordering) 338 of Neumann boundaries for the global problem. 339 340 Not collective 341 342 Input Parameters: 343 + pc - the preconditioning context 344 - NeumannBoundaries - sequential index set defining the subdomain part of Neumann boundaries (can be NULL) 345 346 Level: intermediate 347 348 Notes: 349 350 .seealso: PCBDDC 351 @*/ 352 PetscErrorCode PCBDDCSetNeumannBoundaries(PC pc,IS NeumannBoundaries) 353 { 354 PetscErrorCode ierr; 355 356 PetscFunctionBegin; 357 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 358 PetscValidHeaderSpecific(NeumannBoundaries,IS_CLASSID,2); 359 ierr = PetscTryMethod(pc,"PCBDDCSetNeumannBoundaries_C",(PC,IS),(pc,NeumannBoundaries));CHKERRQ(ierr); 360 PetscFunctionReturn(0); 361 } 362 /* -------------------------------------------------------------------------- */ 363 364 #undef __FUNCT__ 365 #define __FUNCT__ "PCBDDCGetDirichletBoundaries_BDDC" 366 static PetscErrorCode PCBDDCGetDirichletBoundaries_BDDC(PC pc,IS *DirichletBoundaries) 367 { 368 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 369 370 PetscFunctionBegin; 371 *DirichletBoundaries = pcbddc->DirichletBoundaries; 372 PetscFunctionReturn(0); 373 } 374 375 #undef __FUNCT__ 376 #define __FUNCT__ "PCBDDCGetDirichletBoundaries" 377 /*@ 378 PCBDDCGetDirichletBoundaries - Get index set defining subdomain part (in local ordering) 379 of Dirichlet boundaries for the global problem. 380 381 Not collective 382 383 Input Parameters: 384 + pc - the preconditioning context 385 386 Output Parameters: 387 + DirichletBoundaries - index set defining the subdomain part of Dirichlet boundaries 388 389 Level: intermediate 390 391 Notes: 392 393 .seealso: PCBDDC 394 @*/ 395 PetscErrorCode PCBDDCGetDirichletBoundaries(PC pc,IS *DirichletBoundaries) 396 { 397 PetscErrorCode ierr; 398 399 PetscFunctionBegin; 400 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 401 ierr = PetscUseMethod(pc,"PCBDDCGetDirichletBoundaries_C",(PC,IS*),(pc,DirichletBoundaries));CHKERRQ(ierr); 402 PetscFunctionReturn(0); 403 } 404 /* -------------------------------------------------------------------------- */ 405 406 #undef __FUNCT__ 407 #define __FUNCT__ "PCBDDCGetNeumannBoundaries_BDDC" 408 static PetscErrorCode PCBDDCGetNeumannBoundaries_BDDC(PC pc,IS *NeumannBoundaries) 409 { 410 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 411 412 PetscFunctionBegin; 413 *NeumannBoundaries = pcbddc->NeumannBoundaries; 414 PetscFunctionReturn(0); 415 } 416 417 #undef __FUNCT__ 418 #define __FUNCT__ "PCBDDCGetNeumannBoundaries" 419 /*@ 420 PCBDDCGetNeumannBoundaries - Get index set defining subdomain part (in local ordering) 421 of Neumann boundaries for the global problem. 422 423 Not collective 424 425 Input Parameters: 426 + pc - the preconditioning context 427 428 Output Parameters: 429 + NeumannBoundaries - index set defining the subdomain part of Neumann boundaries 430 431 Level: intermediate 432 433 Notes: 434 435 .seealso: PCBDDC 436 @*/ 437 PetscErrorCode PCBDDCGetNeumannBoundaries(PC pc,IS *NeumannBoundaries) 438 { 439 PetscErrorCode ierr; 440 441 PetscFunctionBegin; 442 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 443 ierr = PetscUseMethod(pc,"PCBDDCGetNeumannBoundaries_C",(PC,IS*),(pc,NeumannBoundaries));CHKERRQ(ierr); 444 PetscFunctionReturn(0); 445 } 446 /* -------------------------------------------------------------------------- */ 447 448 #undef __FUNCT__ 449 #define __FUNCT__ "PCBDDCSetLocalAdjacencyGraph_BDDC" 450 static PetscErrorCode PCBDDCSetLocalAdjacencyGraph_BDDC(PC pc, PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode) 451 { 452 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 453 PCBDDCGraph mat_graph = pcbddc->mat_graph; 454 PetscErrorCode ierr; 455 456 PetscFunctionBegin; 457 /* free old CSR */ 458 ierr = PCBDDCGraphResetCSR(mat_graph);CHKERRQ(ierr); 459 /* TODO: PCBDDCGraphSetAdjacency */ 460 /* get CSR into graph structure */ 461 if (copymode == PETSC_COPY_VALUES) { 462 ierr = PetscMalloc((nvtxs+1)*sizeof(PetscInt),&mat_graph->xadj);CHKERRQ(ierr); 463 ierr = PetscMalloc(xadj[nvtxs]*sizeof(PetscInt),&mat_graph->adjncy);CHKERRQ(ierr); 464 ierr = PetscMemcpy(mat_graph->xadj,xadj,(nvtxs+1)*sizeof(PetscInt));CHKERRQ(ierr); 465 ierr = PetscMemcpy(mat_graph->adjncy,adjncy,xadj[nvtxs]*sizeof(PetscInt));CHKERRQ(ierr); 466 } else if (copymode == PETSC_OWN_POINTER) { 467 mat_graph->xadj = (PetscInt*)xadj; 468 mat_graph->adjncy = (PetscInt*)adjncy; 469 } 470 mat_graph->nvtxs_csr = nvtxs; 471 PetscFunctionReturn(0); 472 } 473 474 #undef __FUNCT__ 475 #define __FUNCT__ "PCBDDCSetLocalAdjacencyGraph" 476 /*@ 477 PCBDDCSetLocalAdjacencyGraph - Set CSR graph of local matrix for use of PCBDDC. 478 479 Not collective 480 481 Input Parameters: 482 + pc - the preconditioning context 483 - nvtxs - number of local vertices of the graph 484 - xadj, adjncy - the CSR graph 485 - copymode - either PETSC_COPY_VALUES or PETSC_OWN_POINTER. In the former case the user must free the array passed in; 486 in the latter case, memory must be obtained with PetscMalloc. 487 488 Level: intermediate 489 490 Notes: 491 492 .seealso: PCBDDC 493 @*/ 494 PetscErrorCode PCBDDCSetLocalAdjacencyGraph(PC pc,PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode) 495 { 496 void (*f)(void) = 0; 497 PetscErrorCode ierr; 498 499 PetscFunctionBegin; 500 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 501 PetscValidIntPointer(xadj,3); 502 PetscValidIntPointer(xadj,4); 503 if (copymode != PETSC_COPY_VALUES && copymode != PETSC_OWN_POINTER) { 504 SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported copy mode %d in %s\n",copymode,__FUNCT__); 505 } 506 ierr = PetscTryMethod(pc,"PCBDDCSetLocalAdjacencyGraph_C",(PC,PetscInt,const PetscInt[],const PetscInt[],PetscCopyMode),(pc,nvtxs,xadj,adjncy,copymode));CHKERRQ(ierr); 507 /* free arrays if PCBDDC is not the PC type */ 508 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",&f);CHKERRQ(ierr); 509 if (!f && copymode == PETSC_OWN_POINTER) { 510 ierr = PetscFree(xadj);CHKERRQ(ierr); 511 ierr = PetscFree(adjncy);CHKERRQ(ierr); 512 } 513 PetscFunctionReturn(0); 514 } 515 /* -------------------------------------------------------------------------- */ 516 517 #undef __FUNCT__ 518 #define __FUNCT__ "PCBDDCSetDofsSplitting_BDDC" 519 static PetscErrorCode PCBDDCSetDofsSplitting_BDDC(PC pc,PetscInt n_is, IS ISForDofs[]) 520 { 521 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 522 PetscInt i; 523 PetscErrorCode ierr; 524 525 PetscFunctionBegin; 526 /* Destroy ISes if they were already set */ 527 for (i=0;i<pcbddc->n_ISForDofs;i++) { 528 ierr = ISDestroy(&pcbddc->ISForDofs[i]);CHKERRQ(ierr); 529 } 530 ierr = PetscFree(pcbddc->ISForDofs);CHKERRQ(ierr); 531 /* allocate space then set */ 532 ierr = PetscMalloc(n_is*sizeof(IS),&pcbddc->ISForDofs);CHKERRQ(ierr); 533 for (i=0;i<n_is;i++) { 534 ierr = PetscObjectReference((PetscObject)ISForDofs[i]);CHKERRQ(ierr); 535 pcbddc->ISForDofs[i]=ISForDofs[i]; 536 } 537 pcbddc->n_ISForDofs=n_is; 538 PetscFunctionReturn(0); 539 } 540 541 #undef __FUNCT__ 542 #define __FUNCT__ "PCBDDCSetDofsSplitting" 543 /*@ 544 PCBDDCSetDofsSplitting - Set index sets defining fields of local mat. 545 546 Not collective 547 548 Input Parameters: 549 + pc - the preconditioning context 550 - n - number of index sets defining the fields 551 - IS[] - array of IS describing the fields 552 553 Level: intermediate 554 555 Notes: 556 557 .seealso: PCBDDC 558 @*/ 559 PetscErrorCode PCBDDCSetDofsSplitting(PC pc,PetscInt n_is, IS ISForDofs[]) 560 { 561 PetscInt i; 562 PetscErrorCode ierr; 563 564 PetscFunctionBegin; 565 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 566 for (i=0;i<n_is;i++) { 567 PetscValidHeaderSpecific(ISForDofs[i],IS_CLASSID,2); 568 } 569 ierr = PetscTryMethod(pc,"PCBDDCSetDofsSplitting_C",(PC,PetscInt,IS[]),(pc,n_is,ISForDofs));CHKERRQ(ierr); 570 PetscFunctionReturn(0); 571 } 572 /* -------------------------------------------------------------------------- */ 573 #undef __FUNCT__ 574 #define __FUNCT__ "PCPreSolve_BDDC" 575 /* -------------------------------------------------------------------------- */ 576 /* 577 PCPreSolve_BDDC - Changes the right hand side and (if necessary) the initial 578 guess if a transformation of basis approach has been selected. 579 580 Input Parameter: 581 + pc - the preconditioner contex 582 583 Application Interface Routine: PCPreSolve() 584 585 Notes: 586 The interface routine PCPreSolve() is not usually called directly by 587 the user, but instead is called by KSPSolve(). 588 */ 589 static PetscErrorCode PCPreSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x) 590 { 591 PetscErrorCode ierr; 592 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 593 PC_IS *pcis = (PC_IS*)(pc->data); 594 Mat_IS *matis = (Mat_IS*)pc->pmat->data; 595 Mat temp_mat; 596 IS dirIS; 597 PetscInt dirsize,i,*is_indices; 598 PetscScalar *array_x,*array_diagonal; 599 Vec used_vec; 600 PetscBool guess_nonzero,flg,bddc_has_dirichlet_boundaries; 601 602 PetscFunctionBegin; 603 /* if we are working with cg, one dirichlet solve can be avoided during Krylov iterations */ 604 if (ksp) { 605 PetscBool iscg; 606 ierr = PetscObjectTypeCompare((PetscObject)ksp,KSPCG,&iscg);CHKERRQ(ierr); 607 if (!iscg) { 608 ierr = PCBDDCSetUseExactDirichlet(pc,PETSC_FALSE);CHKERRQ(ierr); 609 } 610 } 611 /* Creates parallel work vectors used in presolve */ 612 if (!pcbddc->original_rhs) { 613 ierr = VecDuplicate(pcis->vec1_global,&pcbddc->original_rhs);CHKERRQ(ierr); 614 } 615 if (!pcbddc->temp_solution) { 616 ierr = VecDuplicate(pcis->vec1_global,&pcbddc->temp_solution);CHKERRQ(ierr); 617 } 618 if (x) { 619 ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr); 620 used_vec = x; 621 } else { 622 ierr = PetscObjectReference((PetscObject)pcbddc->temp_solution);CHKERRQ(ierr); 623 used_vec = pcbddc->temp_solution; 624 ierr = VecSet(used_vec,0.0);CHKERRQ(ierr); 625 } 626 /* hack into ksp data structure PCPreSolve comes earlier in src/ksp/ksp/interface/itfunc.c */ 627 if (ksp) { 628 ierr = KSPGetInitialGuessNonzero(ksp,&guess_nonzero);CHKERRQ(ierr); 629 if (!guess_nonzero) { 630 ierr = VecSet(used_vec,0.0);CHKERRQ(ierr); 631 } 632 } 633 634 /* TODO: remove when Dirichlet boundaries will be shared */ 635 ierr = PCBDDCGetDirichletBoundaries(pc,&dirIS);CHKERRQ(ierr); 636 flg = PETSC_FALSE; 637 if (dirIS) flg = PETSC_TRUE; 638 ierr = MPI_Allreduce(&flg,&bddc_has_dirichlet_boundaries,1,MPIU_BOOL,MPI_LOR,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr); 639 640 /* store the original rhs */ 641 ierr = VecCopy(rhs,pcbddc->original_rhs);CHKERRQ(ierr); 642 643 /* Take into account zeroed rows -> change rhs and store solution removed */ 644 if (rhs && bddc_has_dirichlet_boundaries) { 645 ierr = MatGetDiagonal(pc->pmat,pcis->vec1_global);CHKERRQ(ierr); 646 ierr = VecPointwiseDivide(pcis->vec1_global,rhs,pcis->vec1_global);CHKERRQ(ierr); 647 ierr = VecScatterBegin(matis->ctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 648 ierr = VecScatterEnd(matis->ctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 649 ierr = VecScatterBegin(matis->ctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 650 ierr = VecScatterEnd(matis->ctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 651 if (dirIS) { 652 ierr = ISGetSize(dirIS,&dirsize);CHKERRQ(ierr); 653 ierr = VecGetArray(pcis->vec1_N,&array_x);CHKERRQ(ierr); 654 ierr = VecGetArray(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr); 655 ierr = ISGetIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr); 656 for (i=0; i<dirsize; i++) array_x[is_indices[i]] = array_diagonal[is_indices[i]]; 657 ierr = ISRestoreIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr); 658 ierr = VecRestoreArray(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr); 659 ierr = VecRestoreArray(pcis->vec1_N,&array_x);CHKERRQ(ierr); 660 } 661 ierr = VecScatterBegin(matis->ctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 662 ierr = VecScatterEnd(matis->ctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 663 664 /* remove the computed solution from the rhs */ 665 ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr); 666 ierr = MatMultAdd(pc->pmat,used_vec,rhs,rhs);CHKERRQ(ierr); 667 ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr); 668 } 669 670 /* store partially computed solution and set initial guess */ 671 if (x) { 672 ierr = VecCopy(used_vec,pcbddc->temp_solution);CHKERRQ(ierr); 673 ierr = VecSet(used_vec,0.0);CHKERRQ(ierr); 674 if (pcbddc->use_exact_dirichlet_trick) { 675 ierr = VecScatterBegin(pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 676 ierr = VecScatterEnd (pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 677 ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr); 678 ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 679 ierr = VecScatterEnd (pcis->global_to_D,pcis->vec2_D,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 680 if (ksp) { 681 ierr = KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);CHKERRQ(ierr); 682 } 683 } 684 } 685 686 /* prepare MatMult and rhs for solver */ 687 if (pcbddc->use_change_of_basis) { 688 /* swap pointers for local matrices */ 689 temp_mat = matis->A; 690 matis->A = pcbddc->local_mat; 691 pcbddc->local_mat = temp_mat; 692 if (rhs) { 693 /* Get local rhs and apply transformation of basis */ 694 ierr = VecScatterBegin(pcis->global_to_B,rhs,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 695 ierr = VecScatterEnd (pcis->global_to_B,rhs,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 696 /* from original basis to modified basis */ 697 ierr = MatMultTranspose(pcbddc->ChangeOfBasisMatrix,pcis->vec1_B,pcis->vec2_B);CHKERRQ(ierr); 698 /* put back modified values into the global vec using INSERT_VALUES copy mode */ 699 ierr = VecScatterBegin(pcis->global_to_B,pcis->vec2_B,rhs,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 700 ierr = VecScatterEnd (pcis->global_to_B,pcis->vec2_B,rhs,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 701 } 702 } 703 704 /* remove nullspace if present */ 705 if (ksp && pcbddc->NullSpace) { 706 ierr = MatNullSpaceRemove(pcbddc->NullSpace,used_vec);CHKERRQ(ierr); 707 ierr = MatNullSpaceRemove(pcbddc->NullSpace,rhs);CHKERRQ(ierr); 708 } 709 ierr = VecDestroy(&used_vec);CHKERRQ(ierr); 710 PetscFunctionReturn(0); 711 } 712 /* -------------------------------------------------------------------------- */ 713 #undef __FUNCT__ 714 #define __FUNCT__ "PCPostSolve_BDDC" 715 /* -------------------------------------------------------------------------- */ 716 /* 717 PCPostSolve_BDDC - Changes the computed solution if a transformation of basis 718 approach has been selected. Also, restores rhs to its original state. 719 720 Input Parameter: 721 + pc - the preconditioner contex 722 723 Application Interface Routine: PCPostSolve() 724 725 Notes: 726 The interface routine PCPostSolve() is not usually called directly by 727 the user, but instead is called by KSPSolve(). 728 */ 729 static PetscErrorCode PCPostSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x) 730 { 731 PetscErrorCode ierr; 732 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 733 PC_IS *pcis = (PC_IS*)(pc->data); 734 Mat_IS *matis = (Mat_IS*)pc->pmat->data; 735 Mat temp_mat; 736 737 PetscFunctionBegin; 738 if (pcbddc->use_change_of_basis) { 739 /* swap pointers for local matrices */ 740 temp_mat = matis->A; 741 matis->A = pcbddc->local_mat; 742 pcbddc->local_mat = temp_mat; 743 } 744 if (pcbddc->use_change_of_basis && x) { 745 /* Get Local boundary and apply transformation of basis to solution vector */ 746 ierr = VecScatterBegin(pcis->global_to_B,x,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 747 ierr = VecScatterEnd (pcis->global_to_B,x,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 748 /* from modified basis to original basis */ 749 ierr = MatMult(pcbddc->ChangeOfBasisMatrix,pcis->vec1_B,pcis->vec2_B);CHKERRQ(ierr); 750 /* put back modified values into the global vec using INSERT_VALUES copy mode */ 751 ierr = VecScatterBegin(pcis->global_to_B,pcis->vec2_B,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 752 ierr = VecScatterEnd (pcis->global_to_B,pcis->vec2_B,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 753 } 754 /* add solution removed in presolve */ 755 if (x) { 756 ierr = VecAXPY(x,1.0,pcbddc->temp_solution);CHKERRQ(ierr); 757 } 758 /* restore rhs to its original state */ 759 if (rhs) { 760 ierr = VecCopy(pcbddc->original_rhs,rhs);CHKERRQ(ierr); 761 } 762 PetscFunctionReturn(0); 763 } 764 /* -------------------------------------------------------------------------- */ 765 #undef __FUNCT__ 766 #define __FUNCT__ "PCSetUp_BDDC" 767 /* -------------------------------------------------------------------------- */ 768 /* 769 PCSetUp_BDDC - Prepares for the use of the BDDC preconditioner 770 by setting data structures and options. 771 772 Input Parameter: 773 + pc - the preconditioner context 774 775 Application Interface Routine: PCSetUp() 776 777 Notes: 778 The interface routine PCSetUp() is not usually called directly by 779 the user, but instead is called by PCApply() if necessary. 780 */ 781 PetscErrorCode PCSetUp_BDDC(PC pc) 782 { 783 PetscErrorCode ierr; 784 PC_BDDC* pcbddc = (PC_BDDC*)pc->data; 785 MatStructure flag; 786 PetscBool computeis,computetopography,computesolvers; 787 788 PetscFunctionBegin; 789 /* the following lines of code should be replaced by a better logic between PCIS, PCNN, PCBDDC and other nonoverlapping preconditioners */ 790 /* PCIS does not support MatStructures different from SAME_PRECONDITIONER */ 791 /* For BDDC we need to define a local "Neumann" problem different to that defined in PCISSetup 792 Also, BDDC directly build the Dirichlet problem */ 793 /* Get stdout for dbg */ 794 if (pcbddc->dbg_flag && !pcbddc->dbg_viewer) { 795 ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)pc),&pcbddc->dbg_viewer);CHKERRQ(ierr); 796 ierr = PetscViewerASCIISynchronizedAllow(pcbddc->dbg_viewer,PETSC_TRUE);CHKERRQ(ierr); 797 if (pcbddc->current_level) { 798 ierr = PetscViewerASCIIAddTab(pcbddc->dbg_viewer,2);CHKERRQ(ierr); 799 } 800 } 801 /* first attempt to split work */ 802 if (pc->setupcalled) { 803 computeis = PETSC_FALSE; 804 ierr = PCGetOperators(pc,NULL,NULL,&flag);CHKERRQ(ierr); 805 if (flag == SAME_PRECONDITIONER) { 806 computetopography = PETSC_FALSE; 807 computesolvers = PETSC_FALSE; 808 } else if (flag == SAME_NONZERO_PATTERN) { 809 computetopography = PETSC_FALSE; 810 computesolvers = PETSC_TRUE; 811 } else { /* DIFFERENT_NONZERO_PATTERN */ 812 computetopography = PETSC_TRUE; 813 computesolvers = PETSC_TRUE; 814 } 815 } else { 816 computeis = PETSC_TRUE; 817 computetopography = PETSC_TRUE; 818 computesolvers = PETSC_TRUE; 819 } 820 /* Set up all the "iterative substructuring" common block without computing solvers */ 821 if (computeis) { 822 /* HACK INTO PCIS */ 823 PC_IS* pcis = (PC_IS*)pc->data; 824 pcis->computesolvers = PETSC_FALSE; 825 ierr = PCISSetUp(pc);CHKERRQ(ierr); 826 } 827 /* Analyze interface and set up local constraint and change of basis matrices */ 828 if (computetopography) { 829 /* reset data */ 830 ierr = PCBDDCResetTopography(pc);CHKERRQ(ierr); 831 ierr = PCBDDCAnalyzeInterface(pc);CHKERRQ(ierr); 832 ierr = PCBDDCConstraintsSetUp(pc);CHKERRQ(ierr); 833 } 834 if (computesolvers) { 835 /* reset data */ 836 ierr = PCBDDCResetSolvers(pc);CHKERRQ(ierr); 837 ierr = PCBDDCScalingDestroy(pc);CHKERRQ(ierr); 838 /* Create coarse and local stuffs */ 839 ierr = PCBDDCSetUpSolvers(pc);CHKERRQ(ierr); 840 ierr = PCBDDCScalingSetUp(pc);CHKERRQ(ierr); 841 } 842 if (pcbddc->dbg_flag && pcbddc->current_level) { 843 ierr = PetscViewerASCIISubtractTab(pcbddc->dbg_viewer,2);CHKERRQ(ierr); 844 } 845 PetscFunctionReturn(0); 846 } 847 848 /* -------------------------------------------------------------------------- */ 849 /* 850 PCApply_BDDC - Applies the BDDC preconditioner to a vector. 851 852 Input Parameters: 853 . pc - the preconditioner context 854 . r - input vector (global) 855 856 Output Parameter: 857 . z - output vector (global) 858 859 Application Interface Routine: PCApply() 860 */ 861 #undef __FUNCT__ 862 #define __FUNCT__ "PCApply_BDDC" 863 PetscErrorCode PCApply_BDDC(PC pc,Vec r,Vec z) 864 { 865 PC_IS *pcis = (PC_IS*)(pc->data); 866 PC_BDDC *pcbddc = (PC_BDDC*)(pc->data); 867 PetscErrorCode ierr; 868 const PetscScalar one = 1.0; 869 const PetscScalar m_one = -1.0; 870 const PetscScalar zero = 0.0; 871 872 /* This code is similar to that provided in nn.c for PCNN 873 NN interface preconditioner changed to BDDC 874 Added support for M_3 preconditioner in the reference article (code is active if pcbddc->switch_static = PETSC_TRUE) */ 875 876 PetscFunctionBegin; 877 if (!pcbddc->use_exact_dirichlet_trick) { 878 /* First Dirichlet solve */ 879 ierr = VecScatterBegin(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 880 ierr = VecScatterEnd (pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 881 ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr); 882 /* 883 Assembling right hand side for BDDC operator 884 - pcis->vec1_D for the Dirichlet part (if needed, i.e. prec_flag=PETSC_TRUE) 885 - pcis->vec1_B the interface part of the global vector z 886 */ 887 ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr); 888 ierr = MatMult(pcis->A_BI,pcis->vec2_D,pcis->vec1_B);CHKERRQ(ierr); 889 if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec2_D,pcis->vec1_D,pcis->vec1_D);CHKERRQ(ierr); } 890 ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr); 891 ierr = VecCopy(r,z);CHKERRQ(ierr); 892 ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 893 ierr = VecScatterEnd (pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 894 ierr = PCBDDCScalingRestriction(pc,z,pcis->vec1_B);CHKERRQ(ierr); 895 } else { 896 ierr = VecSet(pcis->vec1_D,zero);CHKERRQ(ierr); 897 ierr = VecSet(pcis->vec2_D,zero);CHKERRQ(ierr); 898 ierr = PCBDDCScalingRestriction(pc,r,pcis->vec1_B);CHKERRQ(ierr); 899 } 900 901 /* Apply interface preconditioner 902 input/output vecs: pcis->vec1_B and pcis->vec1_D */ 903 ierr = PCBDDCApplyInterfacePreconditioner(pc);CHKERRQ(ierr); 904 905 /* Apply transpose of partition of unity operator */ 906 ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,z);CHKERRQ(ierr); 907 908 /* Second Dirichlet solve and assembling of output */ 909 ierr = VecScatterBegin(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 910 ierr = VecScatterEnd (pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 911 ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec3_D);CHKERRQ(ierr); 912 if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec1_D,pcis->vec3_D,pcis->vec3_D);CHKERRQ(ierr); } 913 ierr = KSPSolve(pcbddc->ksp_D,pcis->vec3_D,pcbddc->vec4_D);CHKERRQ(ierr); 914 ierr = VecScale(pcbddc->vec4_D,m_one);CHKERRQ(ierr); 915 if (pcbddc->switch_static) { ierr = VecAXPY (pcbddc->vec4_D,one,pcis->vec1_D);CHKERRQ(ierr); } 916 ierr = VecAXPY (pcis->vec2_D,one,pcbddc->vec4_D);CHKERRQ(ierr); 917 ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 918 ierr = VecScatterEnd (pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 919 PetscFunctionReturn(0); 920 } 921 /* -------------------------------------------------------------------------- */ 922 923 #undef __FUNCT__ 924 #define __FUNCT__ "PCDestroy_BDDC" 925 PetscErrorCode PCDestroy_BDDC(PC pc) 926 { 927 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 928 PetscErrorCode ierr; 929 930 PetscFunctionBegin; 931 /* free data created by PCIS */ 932 ierr = PCISDestroy(pc);CHKERRQ(ierr); 933 /* free BDDC custom data */ 934 ierr = PCBDDCResetCustomization(pc);CHKERRQ(ierr); 935 /* destroy objects related to topography */ 936 ierr = PCBDDCResetTopography(pc);CHKERRQ(ierr); 937 /* free allocated graph structure */ 938 ierr = PetscFree(pcbddc->mat_graph);CHKERRQ(ierr); 939 /* free data for scaling operator */ 940 ierr = PCBDDCScalingDestroy(pc);CHKERRQ(ierr); 941 /* free solvers stuff */ 942 ierr = PCBDDCResetSolvers(pc);CHKERRQ(ierr); 943 ierr = KSPDestroy(&pcbddc->ksp_D);CHKERRQ(ierr); 944 ierr = KSPDestroy(&pcbddc->ksp_R);CHKERRQ(ierr); 945 ierr = KSPDestroy(&pcbddc->coarse_ksp);CHKERRQ(ierr); 946 ierr = MatDestroy(&pcbddc->local_mat);CHKERRQ(ierr); 947 /* free global vectors needed in presolve */ 948 ierr = VecDestroy(&pcbddc->temp_solution);CHKERRQ(ierr); 949 ierr = VecDestroy(&pcbddc->original_rhs);CHKERRQ(ierr); 950 /* remove functions */ 951 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",NULL);CHKERRQ(ierr); 952 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",NULL);CHKERRQ(ierr); 953 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",NULL);CHKERRQ(ierr); 954 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",NULL);CHKERRQ(ierr); 955 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",NULL);CHKERRQ(ierr); 956 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",NULL);CHKERRQ(ierr); 957 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",NULL);CHKERRQ(ierr); 958 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",NULL);CHKERRQ(ierr); 959 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",NULL);CHKERRQ(ierr); 960 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",NULL);CHKERRQ(ierr); 961 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",NULL);CHKERRQ(ierr); 962 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",NULL);CHKERRQ(ierr); 963 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",NULL);CHKERRQ(ierr); 964 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",NULL);CHKERRQ(ierr); 965 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",NULL);CHKERRQ(ierr); 966 /* Free the private data structure */ 967 ierr = PetscFree(pc->data);CHKERRQ(ierr); 968 PetscFunctionReturn(0); 969 } 970 /* -------------------------------------------------------------------------- */ 971 972 #undef __FUNCT__ 973 #define __FUNCT__ "PCBDDCMatFETIDPGetRHS_BDDC" 974 static PetscErrorCode PCBDDCMatFETIDPGetRHS_BDDC(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs) 975 { 976 FETIDPMat_ctx mat_ctx; 977 PC_IS* pcis; 978 PC_BDDC* pcbddc; 979 PetscErrorCode ierr; 980 981 PetscFunctionBegin; 982 ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr); 983 pcis = (PC_IS*)mat_ctx->pc->data; 984 pcbddc = (PC_BDDC*)mat_ctx->pc->data; 985 986 /* change of basis for physical rhs if needed 987 It also changes the rhs in case of dirichlet boundaries */ 988 ierr = PCPreSolve_BDDC(mat_ctx->pc,NULL,standard_rhs,NULL);CHKERRQ(ierr); 989 /* store vectors for computation of fetidp final solution */ 990 ierr = VecScatterBegin(pcis->global_to_D,standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 991 ierr = VecScatterEnd(pcis->global_to_D,standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 992 /* scale rhs since it should be unassembled */ 993 /* TODO use counter scaling? (also below) */ 994 ierr = VecScatterBegin(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 995 ierr = VecScatterEnd(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 996 /* Apply partition of unity */ 997 ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr); 998 /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */ 999 if (!pcbddc->switch_static) { 1000 /* compute partially subassembled Schur complement right-hand side */ 1001 ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr); 1002 ierr = MatMult(pcis->A_BI,pcis->vec1_D,pcis->vec1_B);CHKERRQ(ierr); 1003 ierr = VecAXPY(mat_ctx->temp_solution_B,-1.0,pcis->vec1_B);CHKERRQ(ierr); 1004 ierr = VecSet(standard_rhs,0.0);CHKERRQ(ierr); 1005 ierr = VecScatterBegin(pcis->global_to_B,mat_ctx->temp_solution_B,standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1006 ierr = VecScatterEnd(pcis->global_to_B,mat_ctx->temp_solution_B,standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1007 /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */ 1008 ierr = VecScatterBegin(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1009 ierr = VecScatterEnd(pcis->global_to_B,standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1010 ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr); 1011 } 1012 /* BDDC rhs */ 1013 ierr = VecCopy(mat_ctx->temp_solution_B,pcis->vec1_B);CHKERRQ(ierr); 1014 if (pcbddc->switch_static) { 1015 ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr); 1016 } 1017 /* apply BDDC */ 1018 ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc);CHKERRQ(ierr); 1019 /* Application of B_delta and assembling of rhs for fetidp fluxes */ 1020 ierr = VecSet(fetidp_flux_rhs,0.0);CHKERRQ(ierr); 1021 ierr = MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local);CHKERRQ(ierr); 1022 ierr = VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1023 ierr = VecScatterEnd (mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1024 /* restore original rhs */ 1025 ierr = VecCopy(pcbddc->original_rhs,standard_rhs);CHKERRQ(ierr); 1026 PetscFunctionReturn(0); 1027 } 1028 1029 #undef __FUNCT__ 1030 #define __FUNCT__ "PCBDDCMatFETIDPGetRHS" 1031 /*@ 1032 PCBDDCMatFETIDPGetRHS - Get rhs for FETIDP linear system. 1033 1034 Collective 1035 1036 Input Parameters: 1037 + fetidp_mat - the FETIDP mat obtained by a call to PCBDDCCreateFETIDPOperators 1038 + standard_rhs - the rhs of your linear system 1039 1040 Output Parameters: 1041 + fetidp_flux_rhs - the rhs of the FETIDP linear system 1042 1043 Level: developer 1044 1045 Notes: 1046 1047 .seealso: PCBDDC 1048 @*/ 1049 PetscErrorCode PCBDDCMatFETIDPGetRHS(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs) 1050 { 1051 FETIDPMat_ctx mat_ctx; 1052 PetscErrorCode ierr; 1053 1054 PetscFunctionBegin; 1055 ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr); 1056 ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetRHS_C",(Mat,Vec,Vec),(fetidp_mat,standard_rhs,fetidp_flux_rhs));CHKERRQ(ierr); 1057 PetscFunctionReturn(0); 1058 } 1059 /* -------------------------------------------------------------------------- */ 1060 1061 #undef __FUNCT__ 1062 #define __FUNCT__ "PCBDDCMatFETIDPGetSolution_BDDC" 1063 static PetscErrorCode PCBDDCMatFETIDPGetSolution_BDDC(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol) 1064 { 1065 FETIDPMat_ctx mat_ctx; 1066 PC_IS* pcis; 1067 PC_BDDC* pcbddc; 1068 PetscErrorCode ierr; 1069 1070 PetscFunctionBegin; 1071 ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr); 1072 pcis = (PC_IS*)mat_ctx->pc->data; 1073 pcbddc = (PC_BDDC*)mat_ctx->pc->data; 1074 1075 /* apply B_delta^T */ 1076 ierr = VecScatterBegin(mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1077 ierr = VecScatterEnd (mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1078 ierr = MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr); 1079 /* compute rhs for BDDC application */ 1080 ierr = VecAYPX(pcis->vec1_B,-1.0,mat_ctx->temp_solution_B);CHKERRQ(ierr); 1081 if (pcbddc->switch_static) { 1082 ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr); 1083 } 1084 /* apply BDDC */ 1085 ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc);CHKERRQ(ierr); 1086 /* put values into standard global vector */ 1087 ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1088 ierr = VecScatterEnd (pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1089 if (!pcbddc->switch_static) { 1090 /* compute values into the interior if solved for the partially subassembled Schur complement */ 1091 ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec1_D);CHKERRQ(ierr); 1092 ierr = VecAXPY(mat_ctx->temp_solution_D,-1.0,pcis->vec1_D);CHKERRQ(ierr); 1093 ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr); 1094 } 1095 ierr = VecScatterBegin(pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1096 ierr = VecScatterEnd (pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1097 /* final change of basis if needed 1098 Is also sums the dirichlet part removed during RHS assembling */ 1099 ierr = PCPostSolve_BDDC(mat_ctx->pc,NULL,NULL,standard_sol);CHKERRQ(ierr); 1100 PetscFunctionReturn(0); 1101 } 1102 1103 #undef __FUNCT__ 1104 #define __FUNCT__ "PCBDDCMatFETIDPGetSolution" 1105 /*@ 1106 PCBDDCMatFETIDPGetSolution - Get Solution for FETIDP linear system. 1107 1108 Collective 1109 1110 Input Parameters: 1111 + fetidp_mat - the FETIDP mat obtained by a call to PCBDDCCreateFETIDPOperators 1112 + fetidp_flux_sol - the solution of the FETIDP linear system 1113 1114 Output Parameters: 1115 + standard_sol - the solution on the global domain 1116 1117 Level: developer 1118 1119 Notes: 1120 1121 .seealso: PCBDDC 1122 @*/ 1123 PetscErrorCode PCBDDCMatFETIDPGetSolution(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol) 1124 { 1125 FETIDPMat_ctx mat_ctx; 1126 PetscErrorCode ierr; 1127 1128 PetscFunctionBegin; 1129 ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr); 1130 ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetSolution_C",(Mat,Vec,Vec),(fetidp_mat,fetidp_flux_sol,standard_sol));CHKERRQ(ierr); 1131 PetscFunctionReturn(0); 1132 } 1133 /* -------------------------------------------------------------------------- */ 1134 1135 extern PetscErrorCode FETIDPMatMult(Mat,Vec,Vec); 1136 extern PetscErrorCode PCBDDCDestroyFETIDPMat(Mat); 1137 extern PetscErrorCode FETIDPPCApply(PC,Vec,Vec); 1138 extern PetscErrorCode PCBDDCDestroyFETIDPPC(PC); 1139 1140 #undef __FUNCT__ 1141 #define __FUNCT__ "PCBDDCCreateFETIDPOperators_BDDC" 1142 static PetscErrorCode PCBDDCCreateFETIDPOperators_BDDC(PC pc, Mat *fetidp_mat, PC *fetidp_pc) 1143 { 1144 1145 FETIDPMat_ctx fetidpmat_ctx; 1146 Mat newmat; 1147 FETIDPPC_ctx fetidppc_ctx; 1148 PC newpc; 1149 MPI_Comm comm; 1150 PetscErrorCode ierr; 1151 1152 PetscFunctionBegin; 1153 ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 1154 /* FETIDP linear matrix */ 1155 ierr = PCBDDCCreateFETIDPMatContext(pc,&fetidpmat_ctx);CHKERRQ(ierr); 1156 ierr = PCBDDCSetupFETIDPMatContext(fetidpmat_ctx);CHKERRQ(ierr); 1157 ierr = MatCreateShell(comm,PETSC_DECIDE,PETSC_DECIDE,fetidpmat_ctx->n_lambda,fetidpmat_ctx->n_lambda,fetidpmat_ctx,&newmat);CHKERRQ(ierr); 1158 ierr = MatShellSetOperation(newmat,MATOP_MULT,(void (*)(void))FETIDPMatMult);CHKERRQ(ierr); 1159 ierr = MatShellSetOperation(newmat,MATOP_DESTROY,(void (*)(void))PCBDDCDestroyFETIDPMat);CHKERRQ(ierr); 1160 ierr = MatSetUp(newmat);CHKERRQ(ierr); 1161 /* FETIDP preconditioner */ 1162 ierr = PCBDDCCreateFETIDPPCContext(pc,&fetidppc_ctx);CHKERRQ(ierr); 1163 ierr = PCBDDCSetupFETIDPPCContext(newmat,fetidppc_ctx);CHKERRQ(ierr); 1164 ierr = PCCreate(comm,&newpc);CHKERRQ(ierr); 1165 ierr = PCSetType(newpc,PCSHELL);CHKERRQ(ierr); 1166 ierr = PCShellSetContext(newpc,fetidppc_ctx);CHKERRQ(ierr); 1167 ierr = PCShellSetApply(newpc,FETIDPPCApply);CHKERRQ(ierr); 1168 ierr = PCShellSetDestroy(newpc,PCBDDCDestroyFETIDPPC);CHKERRQ(ierr); 1169 ierr = PCSetOperators(newpc,newmat,newmat,SAME_PRECONDITIONER);CHKERRQ(ierr); 1170 ierr = PCSetUp(newpc);CHKERRQ(ierr); 1171 /* return pointers for objects created */ 1172 *fetidp_mat=newmat; 1173 *fetidp_pc=newpc; 1174 PetscFunctionReturn(0); 1175 } 1176 1177 #undef __FUNCT__ 1178 #define __FUNCT__ "PCBDDCCreateFETIDPOperators" 1179 /*@ 1180 PCBDDCCreateFETIDPOperators - Create operators for FETIDP. 1181 1182 Collective 1183 1184 Input Parameters: 1185 + pc - the BDDC preconditioning context (setup must be already called) 1186 1187 Level: developer 1188 1189 Notes: 1190 1191 .seealso: PCBDDC 1192 @*/ 1193 PetscErrorCode PCBDDCCreateFETIDPOperators(PC pc, Mat *fetidp_mat, PC *fetidp_pc) 1194 { 1195 PetscErrorCode ierr; 1196 1197 PetscFunctionBegin; 1198 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1199 if (pc->setupcalled) { 1200 ierr = PetscUseMethod(pc,"PCBDDCCreateFETIDPOperators_C",(PC,Mat*,PC*),(pc,fetidp_mat,fetidp_pc));CHKERRQ(ierr); 1201 } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"You must call PCSetup_BDDC() first \n"); 1202 PetscFunctionReturn(0); 1203 } 1204 /* -------------------------------------------------------------------------- */ 1205 /*MC 1206 PCBDDC - Balancing Domain Decomposition by Constraints. 1207 1208 Options Database Keys: 1209 . -pcbddc ??? - 1210 1211 Level: intermediate 1212 1213 Notes: The matrix used with this preconditioner must be of type MATIS 1214 1215 Unlike more 'conventional' interface preconditioners, this iterates over ALL the 1216 degrees of freedom, NOT just those on the interface (this allows the use of approximate solvers 1217 on the subdomains). 1218 1219 Options for the coarse grid preconditioner can be set with - 1220 Options for the Dirichlet subproblem can be set with - 1221 Options for the Neumann subproblem can be set with - 1222 1223 Contributed by Stefano Zampini 1224 1225 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, MATIS 1226 M*/ 1227 1228 #undef __FUNCT__ 1229 #define __FUNCT__ "PCCreate_BDDC" 1230 PETSC_EXTERN PetscErrorCode PCCreate_BDDC(PC pc) 1231 { 1232 PetscErrorCode ierr; 1233 PC_BDDC *pcbddc; 1234 1235 PetscFunctionBegin; 1236 /* Creates the private data structure for this preconditioner and attach it to the PC object. */ 1237 ierr = PetscNewLog(pc,PC_BDDC,&pcbddc);CHKERRQ(ierr); 1238 pc->data = (void*)pcbddc; 1239 1240 /* create PCIS data structure */ 1241 ierr = PCISCreate(pc);CHKERRQ(ierr); 1242 1243 /* BDDC customization */ 1244 pcbddc->use_vertices = PETSC_TRUE; 1245 pcbddc->use_edges = PETSC_TRUE; 1246 pcbddc->use_faces = PETSC_FALSE; 1247 pcbddc->use_change_of_basis = PETSC_FALSE; 1248 pcbddc->use_change_on_faces = PETSC_FALSE; 1249 pcbddc->switch_static = PETSC_FALSE; 1250 pcbddc->use_nnsp_true = PETSC_FALSE; /* not yet exposed */ 1251 pcbddc->dbg_flag = 0; 1252 1253 pcbddc->user_primal_vertices = 0; 1254 pcbddc->NullSpace = 0; 1255 pcbddc->temp_solution = 0; 1256 pcbddc->original_rhs = 0; 1257 pcbddc->local_mat = 0; 1258 pcbddc->ChangeOfBasisMatrix = 0; 1259 pcbddc->coarse_vec = 0; 1260 pcbddc->coarse_rhs = 0; 1261 pcbddc->coarse_ksp = 0; 1262 pcbddc->coarse_phi_B = 0; 1263 pcbddc->coarse_phi_D = 0; 1264 pcbddc->coarse_psi_B = 0; 1265 pcbddc->coarse_psi_D = 0; 1266 pcbddc->vec1_P = 0; 1267 pcbddc->vec1_R = 0; 1268 pcbddc->vec2_R = 0; 1269 pcbddc->local_auxmat1 = 0; 1270 pcbddc->local_auxmat2 = 0; 1271 pcbddc->R_to_B = 0; 1272 pcbddc->R_to_D = 0; 1273 pcbddc->ksp_D = 0; 1274 pcbddc->ksp_R = 0; 1275 pcbddc->NeumannBoundaries = 0; 1276 pcbddc->ISForDofs = 0; 1277 pcbddc->ConstraintMatrix = 0; 1278 pcbddc->use_exact_dirichlet_trick = PETSC_TRUE; 1279 pcbddc->coarse_loc_to_glob = 0; 1280 pcbddc->coarsening_ratio = 8; 1281 pcbddc->current_level = 0; 1282 pcbddc->max_levels = 0; 1283 1284 /* create local graph structure */ 1285 ierr = PCBDDCGraphCreate(&pcbddc->mat_graph);CHKERRQ(ierr); 1286 1287 /* scaling */ 1288 pcbddc->use_deluxe_scaling = PETSC_FALSE; 1289 pcbddc->work_scaling = 0; 1290 1291 /* function pointers */ 1292 pc->ops->apply = PCApply_BDDC; 1293 pc->ops->applytranspose = 0; 1294 pc->ops->setup = PCSetUp_BDDC; 1295 pc->ops->destroy = PCDestroy_BDDC; 1296 pc->ops->setfromoptions = PCSetFromOptions_BDDC; 1297 pc->ops->view = 0; 1298 pc->ops->applyrichardson = 0; 1299 pc->ops->applysymmetricleft = 0; 1300 pc->ops->applysymmetricright = 0; 1301 pc->ops->presolve = PCPreSolve_BDDC; 1302 pc->ops->postsolve = PCPostSolve_BDDC; 1303 1304 /* composing function */ 1305 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",PCBDDCSetPrimalVerticesLocalIS_BDDC);CHKERRQ(ierr); 1306 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",PCBDDCSetCoarseningRatio_BDDC);CHKERRQ(ierr); 1307 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",PCBDDCSetLevel_BDDC);CHKERRQ(ierr); 1308 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",PCBDDCSetUseExactDirichlet_BDDC);CHKERRQ(ierr); 1309 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",PCBDDCSetLevels_BDDC);CHKERRQ(ierr); 1310 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",PCBDDCSetNullSpace_BDDC);CHKERRQ(ierr); 1311 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",PCBDDCSetDirichletBoundaries_BDDC);CHKERRQ(ierr); 1312 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",PCBDDCSetNeumannBoundaries_BDDC);CHKERRQ(ierr); 1313 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",PCBDDCGetDirichletBoundaries_BDDC);CHKERRQ(ierr); 1314 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",PCBDDCGetNeumannBoundaries_BDDC);CHKERRQ(ierr); 1315 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",PCBDDCSetDofsSplitting_BDDC);CHKERRQ(ierr); 1316 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",PCBDDCSetLocalAdjacencyGraph_BDDC);CHKERRQ(ierr); 1317 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",PCBDDCCreateFETIDPOperators_BDDC);CHKERRQ(ierr); 1318 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",PCBDDCMatFETIDPGetRHS_BDDC);CHKERRQ(ierr); 1319 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",PCBDDCMatFETIDPGetSolution_BDDC);CHKERRQ(ierr); 1320 PetscFunctionReturn(0); 1321 } 1322 1323