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