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