1 /* TODOLIST 2 3 Solvers 4 - Add support for cholesky for coarse solver (similar to local solvers) 5 - Propagate ksp prefixes for solvers to mat objects? 6 7 User interface 8 - ** DM attached to pc? 9 10 Debugging output 11 - * Better management of verbosity levels of debugging output 12 13 Extra 14 - *** Is it possible to work with PCBDDCGraph on boundary indices only (less memory consumed)? 15 - BDDC with MG framework? 16 17 FETIDP 18 - Move FETIDP code to its own classes 19 20 MATIS related operations contained in BDDC code 21 - Provide general case for subassembling 22 23 */ 24 25 #include <../src/ksp/pc/impls/bddc/bddc.h> /*I "petscpc.h" I*/ /* includes for fortran wrappers */ 26 #include <../src/ksp/pc/impls/bddc/bddcprivate.h> 27 #include <petscblaslapack.h> 28 29 /* temporarily declare it */ 30 PetscErrorCode PCApply_BDDC(PC,Vec,Vec); 31 32 /* -------------------------------------------------------------------------- */ 33 #undef __FUNCT__ 34 #define __FUNCT__ "PCSetFromOptions_BDDC" 35 PetscErrorCode PCSetFromOptions_BDDC(PetscOptions *PetscOptionsObject,PC pc) 36 { 37 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 38 PetscErrorCode ierr; 39 40 PetscFunctionBegin; 41 ierr = PetscOptionsHead(PetscOptionsObject,"BDDC options");CHKERRQ(ierr); 42 /* Verbose debugging */ 43 ierr = PetscOptionsInt("-pc_bddc_check_level","Verbose output for PCBDDC (intended for debug)","none",pcbddc->dbg_flag,&pcbddc->dbg_flag,NULL);CHKERRQ(ierr); 44 /* Primal space cumstomization */ 45 ierr = PetscOptionsBool("-pc_bddc_use_local_mat_graph","Use or not adjacency graph of local mat for interface analysis","none",pcbddc->use_local_adj,&pcbddc->use_local_adj,NULL);CHKERRQ(ierr); 46 ierr = PetscOptionsBool("-pc_bddc_use_vertices","Use or not corner dofs in coarse space","none",pcbddc->use_vertices,&pcbddc->use_vertices,NULL);CHKERRQ(ierr); 47 ierr = PetscOptionsBool("-pc_bddc_use_edges","Use or not edge constraints in coarse space","none",pcbddc->use_edges,&pcbddc->use_edges,NULL);CHKERRQ(ierr); 48 ierr = PetscOptionsBool("-pc_bddc_use_faces","Use or not face constraints in coarse space","none",pcbddc->use_faces,&pcbddc->use_faces,NULL);CHKERRQ(ierr); 49 ierr = PetscOptionsBool("-pc_bddc_use_true_nnsp","Use near null space attached to the matrix without modifications","none",pcbddc->use_nnsp_true,&pcbddc->use_nnsp_true,NULL);CHKERRQ(ierr); 50 ierr = PetscOptionsBool("-pc_bddc_use_qr_single","Use QR factorization for single constraints on cc (QR is used when multiple constraints are present)","none",pcbddc->use_qr_single,&pcbddc->use_qr_single,NULL);CHKERRQ(ierr); 51 /* Change of basis */ 52 ierr = PetscOptionsBool("-pc_bddc_use_change_of_basis","Use or not internal change of basis on local edge nodes","none",pcbddc->use_change_of_basis,&pcbddc->use_change_of_basis,NULL);CHKERRQ(ierr); 53 ierr = PetscOptionsBool("-pc_bddc_use_change_on_faces","Use or not internal change of basis on local face nodes","none",pcbddc->use_change_on_faces,&pcbddc->use_change_on_faces,NULL);CHKERRQ(ierr); 54 if (!pcbddc->use_change_of_basis) { 55 pcbddc->use_change_on_faces = PETSC_FALSE; 56 } 57 /* Switch between M_2 (default) and M_3 preconditioners (as defined by C. Dohrmann in the ref. article) */ 58 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); 59 ierr = PetscOptionsInt("-pc_bddc_coarse_redistribute","Number of procs where to redistribute coarse problem","none",pcbddc->redistribute_coarse,&pcbddc->redistribute_coarse,NULL);CHKERRQ(ierr); 60 ierr = PetscOptionsInt("-pc_bddc_coarsening_ratio","Set coarsening ratio used in multilevel coarsening","none",pcbddc->coarsening_ratio,&pcbddc->coarsening_ratio,NULL);CHKERRQ(ierr); 61 ierr = PetscOptionsInt("-pc_bddc_levels","Set maximum number of levels for multilevel","none",pcbddc->max_levels,&pcbddc->max_levels,NULL);CHKERRQ(ierr); 62 ierr = PetscOptionsBool("-pc_bddc_use_coarse_estimates","Use estimated eigenvalues for coarse problem","none",pcbddc->use_coarse_estimates,&pcbddc->use_coarse_estimates,NULL);CHKERRQ(ierr); 63 ierr = PetscOptionsBool("-pc_bddc_use_deluxe_scaling","Use deluxe scaling for BDDC","none",pcbddc->use_deluxe_scaling,&pcbddc->use_deluxe_scaling,NULL);CHKERRQ(ierr); 64 ierr = PetscOptionsBool("-pc_bddc_schur_rebuild","Whether or not the interface graph for Schur principal minors has to be rebuilt (i.e. define the interface without any adjacency)","none",pcbddc->sub_schurs_rebuild,&pcbddc->sub_schurs_rebuild,NULL);CHKERRQ(ierr); 65 ierr = PetscOptionsInt("-pc_bddc_schur_layers","Number of dofs' layers for the computation of principal minors (i.e. -1 uses all dofs)","none",pcbddc->sub_schurs_layers,&pcbddc->sub_schurs_layers,NULL);CHKERRQ(ierr); 66 ierr = PetscOptionsBool("-pc_bddc_schur_use_useradj","Whether or not the CSR graph specified by the user should be used for computing successive layers (default is to use adj of local mat)","none",pcbddc->sub_schurs_use_useradj,&pcbddc->sub_schurs_use_useradj,NULL);CHKERRQ(ierr); 67 ierr = PetscOptionsBool("-pc_bddc_deluxe_faster","Faster application of deluxe scaling (requires extra work during setup)","none",pcbddc->faster_deluxe,&pcbddc->faster_deluxe,NULL);CHKERRQ(ierr); 68 ierr = PetscOptionsReal("-pc_bddc_adaptive_threshold","Threshold to be used for adaptive selection of constraints","none",pcbddc->adaptive_threshold,&pcbddc->adaptive_threshold,NULL);CHKERRQ(ierr); 69 ierr = PetscOptionsInt("-pc_bddc_adaptive_nmin","Minimum number of constraints per connected components","none",pcbddc->adaptive_nmin,&pcbddc->adaptive_nmin,NULL);CHKERRQ(ierr); 70 ierr = PetscOptionsInt("-pc_bddc_adaptive_nmax","Maximum number of constraints per connected components","none",pcbddc->adaptive_nmax,&pcbddc->adaptive_nmax,NULL);CHKERRQ(ierr); 71 ierr = PetscOptionsBool("-pc_bddc_symmetric","Symmetric computation of primal basis functions","none",pcbddc->symmetric_primal,&pcbddc->symmetric_primal,NULL);CHKERRQ(ierr); 72 ierr = PetscOptionsInt("-pc_bddc_coarse_adj","Number of processors where to map the coarse adjacency list","none",pcbddc->coarse_adj_red,&pcbddc->coarse_adj_red,NULL);CHKERRQ(ierr); 73 ierr = PetscOptionsBool("-pc_bddc_benign_trick","Apply the benign subspace trick to a class of saddle point problems","none",pcbddc->benign_saddle_point,&pcbddc->benign_saddle_point,NULL);CHKERRQ(ierr); 74 ierr = PetscOptionsTail();CHKERRQ(ierr); 75 PetscFunctionReturn(0); 76 } 77 /* -------------------------------------------------------------------------- */ 78 #undef __FUNCT__ 79 #define __FUNCT__ "PCBDDCSetChangeOfBasisMat_BDDC" 80 static PetscErrorCode PCBDDCSetChangeOfBasisMat_BDDC(PC pc, Mat change) 81 { 82 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 83 PetscErrorCode ierr; 84 85 PetscFunctionBegin; 86 ierr = MatDestroy(&pcbddc->user_ChangeOfBasisMatrix);CHKERRQ(ierr); 87 ierr = PetscObjectReference((PetscObject)change);CHKERRQ(ierr); 88 pcbddc->user_ChangeOfBasisMatrix = change; 89 PetscFunctionReturn(0); 90 } 91 #undef __FUNCT__ 92 #define __FUNCT__ "PCBDDCSetChangeOfBasisMat" 93 /*@ 94 PCBDDCSetChangeOfBasisMat - Set user defined change of basis for dofs 95 96 Collective on PC 97 98 Input Parameters: 99 + pc - the preconditioning context 100 - change - the change of basis matrix 101 102 Level: intermediate 103 104 Notes: 105 106 .seealso: PCBDDC 107 @*/ 108 PetscErrorCode PCBDDCSetChangeOfBasisMat(PC pc, Mat change) 109 { 110 PetscErrorCode ierr; 111 112 PetscFunctionBegin; 113 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 114 PetscValidHeaderSpecific(change,MAT_CLASSID,2); 115 PetscCheckSameComm(pc,1,change,2); 116 if (pc->mat) { 117 PetscInt rows_c,cols_c,rows,cols; 118 ierr = MatGetSize(pc->mat,&rows,&cols);CHKERRQ(ierr); 119 ierr = MatGetSize(change,&rows_c,&cols_c);CHKERRQ(ierr); 120 if (rows_c != rows) { 121 SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Invalid number of rows for change of basis matrix! %d != %d",rows_c,rows); 122 } 123 if (cols_c != cols) { 124 SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Invalid number of columns for change of basis matrix! %d != %d",cols_c,cols); 125 } 126 ierr = MatGetLocalSize(pc->mat,&rows,&cols);CHKERRQ(ierr); 127 ierr = MatGetLocalSize(change,&rows_c,&cols_c);CHKERRQ(ierr); 128 if (rows_c != rows) { 129 SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Invalid number of local rows for change of basis matrix! %d != %d",rows_c,rows); 130 } 131 if (cols_c != cols) { 132 SETERRQ2(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Invalid number of local columns for change of basis matrix! %d != %d",cols_c,cols); 133 } 134 } 135 ierr = PetscTryMethod(pc,"PCBDDCSetChangeOfBasisMat_C",(PC,Mat),(pc,change));CHKERRQ(ierr); 136 PetscFunctionReturn(0); 137 } 138 /* -------------------------------------------------------------------------- */ 139 #undef __FUNCT__ 140 #define __FUNCT__ "PCBDDCSetPrimalVerticesLocalIS_BDDC" 141 static PetscErrorCode PCBDDCSetPrimalVerticesLocalIS_BDDC(PC pc, IS PrimalVertices) 142 { 143 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 144 PetscErrorCode ierr; 145 146 PetscFunctionBegin; 147 ierr = ISDestroy(&pcbddc->user_primal_vertices);CHKERRQ(ierr); 148 ierr = PetscObjectReference((PetscObject)PrimalVertices);CHKERRQ(ierr); 149 pcbddc->user_primal_vertices = PrimalVertices; 150 PetscFunctionReturn(0); 151 } 152 #undef __FUNCT__ 153 #define __FUNCT__ "PCBDDCSetPrimalVerticesLocalIS" 154 /*@ 155 PCBDDCSetPrimalVerticesLocalIS - Set additional user defined primal vertices in PCBDDC 156 157 Collective 158 159 Input Parameters: 160 + pc - the preconditioning context 161 - PrimalVertices - index set of primal vertices in local numbering (can be empty) 162 163 Level: intermediate 164 165 Notes: 166 167 .seealso: PCBDDC 168 @*/ 169 PetscErrorCode PCBDDCSetPrimalVerticesLocalIS(PC pc, IS PrimalVertices) 170 { 171 PetscErrorCode ierr; 172 173 PetscFunctionBegin; 174 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 175 PetscValidHeaderSpecific(PrimalVertices,IS_CLASSID,2); 176 PetscCheckSameComm(pc,1,PrimalVertices,2); 177 ierr = PetscTryMethod(pc,"PCBDDCSetPrimalVerticesLocalIS_C",(PC,IS),(pc,PrimalVertices));CHKERRQ(ierr); 178 PetscFunctionReturn(0); 179 } 180 /* -------------------------------------------------------------------------- */ 181 #undef __FUNCT__ 182 #define __FUNCT__ "PCBDDCSetCoarseningRatio_BDDC" 183 static PetscErrorCode PCBDDCSetCoarseningRatio_BDDC(PC pc,PetscInt k) 184 { 185 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 186 187 PetscFunctionBegin; 188 pcbddc->coarsening_ratio = k; 189 PetscFunctionReturn(0); 190 } 191 192 #undef __FUNCT__ 193 #define __FUNCT__ "PCBDDCSetCoarseningRatio" 194 /*@ 195 PCBDDCSetCoarseningRatio - Set coarsening ratio used in multilevel 196 197 Logically collective on PC 198 199 Input Parameters: 200 + pc - the preconditioning context 201 - k - coarsening ratio (H/h at the coarser level) 202 203 Options Database Keys: 204 . -pc_bddc_coarsening_ratio 205 206 Level: intermediate 207 208 Notes: 209 Approximatively k subdomains at the finer level will be aggregated into a single subdomain at the coarser level 210 211 .seealso: PCBDDC, PCBDDCSetLevels() 212 @*/ 213 PetscErrorCode PCBDDCSetCoarseningRatio(PC pc,PetscInt k) 214 { 215 PetscErrorCode ierr; 216 217 PetscFunctionBegin; 218 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 219 PetscValidLogicalCollectiveInt(pc,k,2); 220 ierr = PetscTryMethod(pc,"PCBDDCSetCoarseningRatio_C",(PC,PetscInt),(pc,k));CHKERRQ(ierr); 221 PetscFunctionReturn(0); 222 } 223 224 /* The following functions (PCBDDCSetUseExactDirichlet PCBDDCSetLevel) are not public */ 225 #undef __FUNCT__ 226 #define __FUNCT__ "PCBDDCSetUseExactDirichlet_BDDC" 227 static PetscErrorCode PCBDDCSetUseExactDirichlet_BDDC(PC pc,PetscBool flg) 228 { 229 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 230 231 PetscFunctionBegin; 232 pcbddc->use_exact_dirichlet_trick = flg; 233 PetscFunctionReturn(0); 234 } 235 236 #undef __FUNCT__ 237 #define __FUNCT__ "PCBDDCSetUseExactDirichlet" 238 PetscErrorCode PCBDDCSetUseExactDirichlet(PC pc,PetscBool flg) 239 { 240 PetscErrorCode ierr; 241 242 PetscFunctionBegin; 243 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 244 PetscValidLogicalCollectiveBool(pc,flg,2); 245 ierr = PetscTryMethod(pc,"PCBDDCSetUseExactDirichlet_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr); 246 PetscFunctionReturn(0); 247 } 248 249 #undef __FUNCT__ 250 #define __FUNCT__ "PCBDDCSetLevel_BDDC" 251 static PetscErrorCode PCBDDCSetLevel_BDDC(PC pc,PetscInt level) 252 { 253 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 254 255 PetscFunctionBegin; 256 pcbddc->current_level = level; 257 PetscFunctionReturn(0); 258 } 259 260 #undef __FUNCT__ 261 #define __FUNCT__ "PCBDDCSetLevel" 262 PetscErrorCode PCBDDCSetLevel(PC pc,PetscInt level) 263 { 264 PetscErrorCode ierr; 265 266 PetscFunctionBegin; 267 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 268 PetscValidLogicalCollectiveInt(pc,level,2); 269 ierr = PetscTryMethod(pc,"PCBDDCSetLevel_C",(PC,PetscInt),(pc,level));CHKERRQ(ierr); 270 PetscFunctionReturn(0); 271 } 272 273 #undef __FUNCT__ 274 #define __FUNCT__ "PCBDDCSetLevels_BDDC" 275 static PetscErrorCode PCBDDCSetLevels_BDDC(PC pc,PetscInt levels) 276 { 277 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 278 279 PetscFunctionBegin; 280 pcbddc->max_levels = levels; 281 PetscFunctionReturn(0); 282 } 283 284 #undef __FUNCT__ 285 #define __FUNCT__ "PCBDDCSetLevels" 286 /*@ 287 PCBDDCSetLevels - Sets the maximum number of levels for multilevel 288 289 Logically collective on PC 290 291 Input Parameters: 292 + pc - the preconditioning context 293 - levels - the maximum number of levels (max 9) 294 295 Options Database Keys: 296 . -pc_bddc_levels 297 298 Level: intermediate 299 300 Notes: 301 Default value is 0, i.e. traditional one-level BDDC 302 303 .seealso: PCBDDC, PCBDDCSetCoarseningRatio() 304 @*/ 305 PetscErrorCode PCBDDCSetLevels(PC pc,PetscInt levels) 306 { 307 PetscErrorCode ierr; 308 309 PetscFunctionBegin; 310 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 311 PetscValidLogicalCollectiveInt(pc,levels,2); 312 if (levels > 99) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Maximum number of levels for bddc is 99\n"); 313 ierr = PetscTryMethod(pc,"PCBDDCSetLevels_C",(PC,PetscInt),(pc,levels));CHKERRQ(ierr); 314 PetscFunctionReturn(0); 315 } 316 /* -------------------------------------------------------------------------- */ 317 318 #undef __FUNCT__ 319 #define __FUNCT__ "PCBDDCSetNullSpace_BDDC" 320 static PetscErrorCode PCBDDCSetNullSpace_BDDC(PC pc,MatNullSpace NullSpace) 321 { 322 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 323 PetscErrorCode ierr; 324 325 PetscFunctionBegin; 326 ierr = PetscObjectReference((PetscObject)NullSpace);CHKERRQ(ierr); 327 ierr = MatNullSpaceDestroy(&pcbddc->NullSpace);CHKERRQ(ierr); 328 pcbddc->NullSpace = NullSpace; 329 PetscFunctionReturn(0); 330 } 331 332 #undef __FUNCT__ 333 #define __FUNCT__ "PCBDDCSetNullSpace" 334 /*@ 335 PCBDDCSetNullSpace - Set nullspace for BDDC operator 336 337 Logically collective on PC and MatNullSpace 338 339 Input Parameters: 340 + pc - the preconditioning context 341 - NullSpace - Null space of the linear operator to be preconditioned (Pmat) 342 343 Level: intermediate 344 345 Notes: 346 347 .seealso: PCBDDC 348 @*/ 349 PetscErrorCode PCBDDCSetNullSpace(PC pc,MatNullSpace NullSpace) 350 { 351 PetscErrorCode ierr; 352 353 PetscFunctionBegin; 354 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 355 PetscValidHeaderSpecific(NullSpace,MAT_NULLSPACE_CLASSID,2); 356 PetscCheckSameComm(pc,1,NullSpace,2); 357 ierr = PetscTryMethod(pc,"PCBDDCSetNullSpace_C",(PC,MatNullSpace),(pc,NullSpace));CHKERRQ(ierr); 358 PetscFunctionReturn(0); 359 } 360 /* -------------------------------------------------------------------------- */ 361 362 #undef __FUNCT__ 363 #define __FUNCT__ "PCBDDCSetDirichletBoundaries_BDDC" 364 static PetscErrorCode PCBDDCSetDirichletBoundaries_BDDC(PC pc,IS DirichletBoundaries) 365 { 366 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 367 PetscErrorCode ierr; 368 369 PetscFunctionBegin; 370 /* last user setting takes precendence -> destroy any other customization */ 371 ierr = ISDestroy(&pcbddc->DirichletBoundariesLocal);CHKERRQ(ierr); 372 ierr = ISDestroy(&pcbddc->DirichletBoundaries);CHKERRQ(ierr); 373 ierr = PetscObjectReference((PetscObject)DirichletBoundaries);CHKERRQ(ierr); 374 pcbddc->DirichletBoundaries = DirichletBoundaries; 375 pcbddc->recompute_topography = PETSC_TRUE; 376 PetscFunctionReturn(0); 377 } 378 379 #undef __FUNCT__ 380 #define __FUNCT__ "PCBDDCSetDirichletBoundaries" 381 /*@ 382 PCBDDCSetDirichletBoundaries - Set IS defining Dirichlet boundaries for the global problem. 383 384 Collective 385 386 Input Parameters: 387 + pc - the preconditioning context 388 - DirichletBoundaries - parallel IS defining the Dirichlet boundaries 389 390 Level: intermediate 391 392 Notes: 393 Provide the information if you used MatZeroRows/Columns routines. Any process can list any global node 394 395 .seealso: PCBDDC, PCBDDCSetDirichletBoundariesLocal(), MatZeroRows(), MatZeroRowsColumns() 396 @*/ 397 PetscErrorCode PCBDDCSetDirichletBoundaries(PC pc,IS DirichletBoundaries) 398 { 399 PetscErrorCode ierr; 400 401 PetscFunctionBegin; 402 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 403 PetscValidHeaderSpecific(DirichletBoundaries,IS_CLASSID,2); 404 PetscCheckSameComm(pc,1,DirichletBoundaries,2); 405 ierr = PetscTryMethod(pc,"PCBDDCSetDirichletBoundaries_C",(PC,IS),(pc,DirichletBoundaries));CHKERRQ(ierr); 406 PetscFunctionReturn(0); 407 } 408 /* -------------------------------------------------------------------------- */ 409 410 #undef __FUNCT__ 411 #define __FUNCT__ "PCBDDCSetDirichletBoundariesLocal_BDDC" 412 static PetscErrorCode PCBDDCSetDirichletBoundariesLocal_BDDC(PC pc,IS DirichletBoundaries) 413 { 414 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 415 PetscErrorCode ierr; 416 417 PetscFunctionBegin; 418 /* last user setting takes precendence -> destroy any other customization */ 419 ierr = ISDestroy(&pcbddc->DirichletBoundariesLocal);CHKERRQ(ierr); 420 ierr = ISDestroy(&pcbddc->DirichletBoundaries);CHKERRQ(ierr); 421 ierr = PetscObjectReference((PetscObject)DirichletBoundaries);CHKERRQ(ierr); 422 pcbddc->DirichletBoundariesLocal = DirichletBoundaries; 423 pcbddc->recompute_topography = PETSC_TRUE; 424 PetscFunctionReturn(0); 425 } 426 427 #undef __FUNCT__ 428 #define __FUNCT__ "PCBDDCSetDirichletBoundariesLocal" 429 /*@ 430 PCBDDCSetDirichletBoundariesLocal - Set IS defining Dirichlet boundaries for the global problem in local ordering. 431 432 Collective 433 434 Input Parameters: 435 + pc - the preconditioning context 436 - DirichletBoundaries - parallel IS defining the Dirichlet boundaries (in local ordering) 437 438 Level: intermediate 439 440 Notes: 441 442 .seealso: PCBDDC, PCBDDCSetDirichletBoundaries(), MatZeroRows(), MatZeroRowsColumns() 443 @*/ 444 PetscErrorCode PCBDDCSetDirichletBoundariesLocal(PC pc,IS DirichletBoundaries) 445 { 446 PetscErrorCode ierr; 447 448 PetscFunctionBegin; 449 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 450 PetscValidHeaderSpecific(DirichletBoundaries,IS_CLASSID,2); 451 PetscCheckSameComm(pc,1,DirichletBoundaries,2); 452 ierr = PetscTryMethod(pc,"PCBDDCSetDirichletBoundariesLocal_C",(PC,IS),(pc,DirichletBoundaries));CHKERRQ(ierr); 453 PetscFunctionReturn(0); 454 } 455 /* -------------------------------------------------------------------------- */ 456 457 #undef __FUNCT__ 458 #define __FUNCT__ "PCBDDCSetNeumannBoundaries_BDDC" 459 static PetscErrorCode PCBDDCSetNeumannBoundaries_BDDC(PC pc,IS NeumannBoundaries) 460 { 461 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 462 PetscErrorCode ierr; 463 464 PetscFunctionBegin; 465 /* last user setting takes precendence -> destroy any other customization */ 466 ierr = ISDestroy(&pcbddc->NeumannBoundariesLocal);CHKERRQ(ierr); 467 ierr = ISDestroy(&pcbddc->NeumannBoundaries);CHKERRQ(ierr); 468 ierr = PetscObjectReference((PetscObject)NeumannBoundaries);CHKERRQ(ierr); 469 pcbddc->NeumannBoundaries = NeumannBoundaries; 470 pcbddc->recompute_topography = PETSC_TRUE; 471 PetscFunctionReturn(0); 472 } 473 474 #undef __FUNCT__ 475 #define __FUNCT__ "PCBDDCSetNeumannBoundaries" 476 /*@ 477 PCBDDCSetNeumannBoundaries - Set IS defining Neumann boundaries for the global problem. 478 479 Collective 480 481 Input Parameters: 482 + pc - the preconditioning context 483 - NeumannBoundaries - parallel IS defining the Neumann boundaries 484 485 Level: intermediate 486 487 Notes: 488 Any process can list any global node 489 490 .seealso: PCBDDC, PCBDDCSetNeumannBoundariesLocal() 491 @*/ 492 PetscErrorCode PCBDDCSetNeumannBoundaries(PC pc,IS NeumannBoundaries) 493 { 494 PetscErrorCode ierr; 495 496 PetscFunctionBegin; 497 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 498 PetscValidHeaderSpecific(NeumannBoundaries,IS_CLASSID,2); 499 PetscCheckSameComm(pc,1,NeumannBoundaries,2); 500 ierr = PetscTryMethod(pc,"PCBDDCSetNeumannBoundaries_C",(PC,IS),(pc,NeumannBoundaries));CHKERRQ(ierr); 501 PetscFunctionReturn(0); 502 } 503 /* -------------------------------------------------------------------------- */ 504 505 #undef __FUNCT__ 506 #define __FUNCT__ "PCBDDCSetNeumannBoundariesLocal_BDDC" 507 static PetscErrorCode PCBDDCSetNeumannBoundariesLocal_BDDC(PC pc,IS NeumannBoundaries) 508 { 509 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 510 PetscErrorCode ierr; 511 512 PetscFunctionBegin; 513 /* last user setting takes precendence -> destroy any other customization */ 514 ierr = ISDestroy(&pcbddc->NeumannBoundariesLocal);CHKERRQ(ierr); 515 ierr = ISDestroy(&pcbddc->NeumannBoundaries);CHKERRQ(ierr); 516 ierr = PetscObjectReference((PetscObject)NeumannBoundaries);CHKERRQ(ierr); 517 pcbddc->NeumannBoundariesLocal = NeumannBoundaries; 518 pcbddc->recompute_topography = PETSC_TRUE; 519 PetscFunctionReturn(0); 520 } 521 522 #undef __FUNCT__ 523 #define __FUNCT__ "PCBDDCSetNeumannBoundariesLocal" 524 /*@ 525 PCBDDCSetNeumannBoundariesLocal - Set IS defining Neumann boundaries for the global problem in local ordering. 526 527 Collective 528 529 Input Parameters: 530 + pc - the preconditioning context 531 - NeumannBoundaries - parallel IS defining the subdomain part of Neumann boundaries (in local ordering) 532 533 Level: intermediate 534 535 Notes: 536 537 .seealso: PCBDDC, PCBDDCSetNeumannBoundaries() 538 @*/ 539 PetscErrorCode PCBDDCSetNeumannBoundariesLocal(PC pc,IS NeumannBoundaries) 540 { 541 PetscErrorCode ierr; 542 543 PetscFunctionBegin; 544 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 545 PetscValidHeaderSpecific(NeumannBoundaries,IS_CLASSID,2); 546 PetscCheckSameComm(pc,1,NeumannBoundaries,2); 547 ierr = PetscTryMethod(pc,"PCBDDCSetNeumannBoundariesLocal_C",(PC,IS),(pc,NeumannBoundaries));CHKERRQ(ierr); 548 PetscFunctionReturn(0); 549 } 550 /* -------------------------------------------------------------------------- */ 551 552 #undef __FUNCT__ 553 #define __FUNCT__ "PCBDDCGetDirichletBoundaries_BDDC" 554 static PetscErrorCode PCBDDCGetDirichletBoundaries_BDDC(PC pc,IS *DirichletBoundaries) 555 { 556 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 557 558 PetscFunctionBegin; 559 *DirichletBoundaries = pcbddc->DirichletBoundaries; 560 PetscFunctionReturn(0); 561 } 562 563 #undef __FUNCT__ 564 #define __FUNCT__ "PCBDDCGetDirichletBoundaries" 565 /*@ 566 PCBDDCGetDirichletBoundaries - Get parallel IS for Dirichlet boundaries 567 568 Collective 569 570 Input Parameters: 571 . pc - the preconditioning context 572 573 Output Parameters: 574 . DirichletBoundaries - index set defining the Dirichlet boundaries 575 576 Level: intermediate 577 578 Notes: 579 The IS returned (if any) is the same passed in earlier by the user with PCBDDCSetDirichletBoundaries 580 581 .seealso: PCBDDC 582 @*/ 583 PetscErrorCode PCBDDCGetDirichletBoundaries(PC pc,IS *DirichletBoundaries) 584 { 585 PetscErrorCode ierr; 586 587 PetscFunctionBegin; 588 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 589 ierr = PetscUseMethod(pc,"PCBDDCGetDirichletBoundaries_C",(PC,IS*),(pc,DirichletBoundaries));CHKERRQ(ierr); 590 PetscFunctionReturn(0); 591 } 592 /* -------------------------------------------------------------------------- */ 593 594 #undef __FUNCT__ 595 #define __FUNCT__ "PCBDDCGetDirichletBoundariesLocal_BDDC" 596 static PetscErrorCode PCBDDCGetDirichletBoundariesLocal_BDDC(PC pc,IS *DirichletBoundaries) 597 { 598 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 599 600 PetscFunctionBegin; 601 *DirichletBoundaries = pcbddc->DirichletBoundariesLocal; 602 PetscFunctionReturn(0); 603 } 604 605 #undef __FUNCT__ 606 #define __FUNCT__ "PCBDDCGetDirichletBoundariesLocal" 607 /*@ 608 PCBDDCGetDirichletBoundariesLocal - Get parallel IS for Dirichlet boundaries (in local ordering) 609 610 Collective 611 612 Input Parameters: 613 . pc - the preconditioning context 614 615 Output Parameters: 616 . DirichletBoundaries - index set defining the subdomain part of Dirichlet boundaries 617 618 Level: intermediate 619 620 Notes: 621 The IS returned could be the same passed in earlier by the user (if provided with PCBDDCSetDirichletBoundariesLocal) or a global-to-local map of the global IS (if provided with PCBDDCSetDirichletBoundaries). 622 In the latter case, the IS will be available after PCSetUp. 623 624 .seealso: PCBDDC 625 @*/ 626 PetscErrorCode PCBDDCGetDirichletBoundariesLocal(PC pc,IS *DirichletBoundaries) 627 { 628 PetscErrorCode ierr; 629 630 PetscFunctionBegin; 631 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 632 ierr = PetscUseMethod(pc,"PCBDDCGetDirichletBoundariesLocal_C",(PC,IS*),(pc,DirichletBoundaries));CHKERRQ(ierr); 633 PetscFunctionReturn(0); 634 } 635 /* -------------------------------------------------------------------------- */ 636 637 #undef __FUNCT__ 638 #define __FUNCT__ "PCBDDCGetNeumannBoundaries_BDDC" 639 static PetscErrorCode PCBDDCGetNeumannBoundaries_BDDC(PC pc,IS *NeumannBoundaries) 640 { 641 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 642 643 PetscFunctionBegin; 644 *NeumannBoundaries = pcbddc->NeumannBoundaries; 645 PetscFunctionReturn(0); 646 } 647 648 #undef __FUNCT__ 649 #define __FUNCT__ "PCBDDCGetNeumannBoundaries" 650 /*@ 651 PCBDDCGetNeumannBoundaries - Get parallel IS for Neumann boundaries 652 653 Collective 654 655 Input Parameters: 656 . pc - the preconditioning context 657 658 Output Parameters: 659 . NeumannBoundaries - index set defining the Neumann boundaries 660 661 Level: intermediate 662 663 Notes: 664 The IS returned (if any) is the same passed in earlier by the user with PCBDDCSetNeumannBoundaries 665 666 .seealso: PCBDDC 667 @*/ 668 PetscErrorCode PCBDDCGetNeumannBoundaries(PC pc,IS *NeumannBoundaries) 669 { 670 PetscErrorCode ierr; 671 672 PetscFunctionBegin; 673 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 674 ierr = PetscUseMethod(pc,"PCBDDCGetNeumannBoundaries_C",(PC,IS*),(pc,NeumannBoundaries));CHKERRQ(ierr); 675 PetscFunctionReturn(0); 676 } 677 /* -------------------------------------------------------------------------- */ 678 679 #undef __FUNCT__ 680 #define __FUNCT__ "PCBDDCGetNeumannBoundariesLocal_BDDC" 681 static PetscErrorCode PCBDDCGetNeumannBoundariesLocal_BDDC(PC pc,IS *NeumannBoundaries) 682 { 683 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 684 685 PetscFunctionBegin; 686 *NeumannBoundaries = pcbddc->NeumannBoundariesLocal; 687 PetscFunctionReturn(0); 688 } 689 690 #undef __FUNCT__ 691 #define __FUNCT__ "PCBDDCGetNeumannBoundariesLocal" 692 /*@ 693 PCBDDCGetNeumannBoundariesLocal - Get parallel IS for Neumann boundaries (in local ordering) 694 695 Collective 696 697 Input Parameters: 698 . pc - the preconditioning context 699 700 Output Parameters: 701 . NeumannBoundaries - index set defining the subdomain part of Neumann boundaries 702 703 Level: intermediate 704 705 Notes: 706 The IS returned could be the same passed in earlier by the user (if provided with PCBDDCSetNeumannBoundariesLocal) or a global-to-local map of the global IS (if provided with PCBDDCSetNeumannBoundaries). 707 In the latter case, the IS will be available after PCSetUp. 708 709 .seealso: PCBDDC 710 @*/ 711 PetscErrorCode PCBDDCGetNeumannBoundariesLocal(PC pc,IS *NeumannBoundaries) 712 { 713 PetscErrorCode ierr; 714 715 PetscFunctionBegin; 716 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 717 ierr = PetscUseMethod(pc,"PCBDDCGetNeumannBoundariesLocal_C",(PC,IS*),(pc,NeumannBoundaries));CHKERRQ(ierr); 718 PetscFunctionReturn(0); 719 } 720 /* -------------------------------------------------------------------------- */ 721 722 #undef __FUNCT__ 723 #define __FUNCT__ "PCBDDCSetLocalAdjacencyGraph_BDDC" 724 static PetscErrorCode PCBDDCSetLocalAdjacencyGraph_BDDC(PC pc, PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode) 725 { 726 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 727 PCBDDCGraph mat_graph = pcbddc->mat_graph; 728 PetscErrorCode ierr; 729 730 PetscFunctionBegin; 731 /* free old CSR */ 732 ierr = PCBDDCGraphResetCSR(mat_graph);CHKERRQ(ierr); 733 /* TODO: PCBDDCGraphSetAdjacency */ 734 /* get CSR into graph structure */ 735 if (copymode == PETSC_COPY_VALUES) { 736 ierr = PetscMalloc1(nvtxs+1,&mat_graph->xadj);CHKERRQ(ierr); 737 ierr = PetscMalloc1(xadj[nvtxs],&mat_graph->adjncy);CHKERRQ(ierr); 738 ierr = PetscMemcpy(mat_graph->xadj,xadj,(nvtxs+1)*sizeof(PetscInt));CHKERRQ(ierr); 739 ierr = PetscMemcpy(mat_graph->adjncy,adjncy,xadj[nvtxs]*sizeof(PetscInt));CHKERRQ(ierr); 740 } else if (copymode == PETSC_OWN_POINTER) { 741 mat_graph->xadj = (PetscInt*)xadj; 742 mat_graph->adjncy = (PetscInt*)adjncy; 743 } 744 mat_graph->nvtxs_csr = nvtxs; 745 PetscFunctionReturn(0); 746 } 747 748 #undef __FUNCT__ 749 #define __FUNCT__ "PCBDDCSetLocalAdjacencyGraph" 750 /*@ 751 PCBDDCSetLocalAdjacencyGraph - Set adjacency structure (CSR graph) of the local matrix 752 753 Not collective 754 755 Input Parameters: 756 + pc - the preconditioning context 757 . nvtxs - number of local vertices of the graph (i.e., the size of the local problem) 758 . xadj, adjncy - the CSR graph 759 - copymode - either PETSC_COPY_VALUES or PETSC_OWN_POINTER. 760 761 Level: intermediate 762 763 Notes: 764 765 .seealso: PCBDDC,PetscCopyMode 766 @*/ 767 PetscErrorCode PCBDDCSetLocalAdjacencyGraph(PC pc,PetscInt nvtxs,const PetscInt xadj[],const PetscInt adjncy[], PetscCopyMode copymode) 768 { 769 void (*f)(void) = 0; 770 PetscErrorCode ierr; 771 772 PetscFunctionBegin; 773 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 774 PetscValidIntPointer(xadj,3); 775 PetscValidIntPointer(adjncy,4); 776 if (copymode != PETSC_COPY_VALUES && copymode != PETSC_OWN_POINTER) { 777 SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported copy mode %d in %s\n",copymode,__FUNCT__); 778 } 779 ierr = PetscTryMethod(pc,"PCBDDCSetLocalAdjacencyGraph_C",(PC,PetscInt,const PetscInt[],const PetscInt[],PetscCopyMode),(pc,nvtxs,xadj,adjncy,copymode));CHKERRQ(ierr); 780 /* free arrays if PCBDDC is not the PC type */ 781 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",&f);CHKERRQ(ierr); 782 if (!f && copymode == PETSC_OWN_POINTER) { 783 ierr = PetscFree(xadj);CHKERRQ(ierr); 784 ierr = PetscFree(adjncy);CHKERRQ(ierr); 785 } 786 PetscFunctionReturn(0); 787 } 788 /* -------------------------------------------------------------------------- */ 789 790 #undef __FUNCT__ 791 #define __FUNCT__ "PCBDDCSetDofsSplittingLocal_BDDC" 792 static PetscErrorCode PCBDDCSetDofsSplittingLocal_BDDC(PC pc,PetscInt n_is, IS ISForDofs[]) 793 { 794 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 795 PetscInt i; 796 PetscErrorCode ierr; 797 798 PetscFunctionBegin; 799 /* Destroy ISes if they were already set */ 800 for (i=0;i<pcbddc->n_ISForDofsLocal;i++) { 801 ierr = ISDestroy(&pcbddc->ISForDofsLocal[i]);CHKERRQ(ierr); 802 } 803 ierr = PetscFree(pcbddc->ISForDofsLocal);CHKERRQ(ierr); 804 /* last user setting takes precendence -> destroy any other customization */ 805 for (i=0;i<pcbddc->n_ISForDofs;i++) { 806 ierr = ISDestroy(&pcbddc->ISForDofs[i]);CHKERRQ(ierr); 807 } 808 ierr = PetscFree(pcbddc->ISForDofs);CHKERRQ(ierr); 809 pcbddc->n_ISForDofs = 0; 810 /* allocate space then set */ 811 if (n_is) { 812 ierr = PetscMalloc1(n_is,&pcbddc->ISForDofsLocal);CHKERRQ(ierr); 813 } 814 for (i=0;i<n_is;i++) { 815 ierr = PetscObjectReference((PetscObject)ISForDofs[i]);CHKERRQ(ierr); 816 pcbddc->ISForDofsLocal[i]=ISForDofs[i]; 817 } 818 pcbddc->n_ISForDofsLocal=n_is; 819 if (n_is) pcbddc->user_provided_isfordofs = PETSC_TRUE; 820 pcbddc->recompute_topography = PETSC_TRUE; 821 PetscFunctionReturn(0); 822 } 823 824 #undef __FUNCT__ 825 #define __FUNCT__ "PCBDDCSetDofsSplittingLocal" 826 /*@ 827 PCBDDCSetDofsSplittingLocal - Set index sets defining fields of the local subdomain matrix 828 829 Collective 830 831 Input Parameters: 832 + pc - the preconditioning context 833 . n_is - number of index sets defining the fields 834 - ISForDofs - array of IS describing the fields in local ordering 835 836 Level: intermediate 837 838 Notes: 839 n_is should be the same among processes. Not all nodes need to be listed: unlisted nodes will belong to the complement field. 840 841 .seealso: PCBDDC 842 @*/ 843 PetscErrorCode PCBDDCSetDofsSplittingLocal(PC pc,PetscInt n_is, IS ISForDofs[]) 844 { 845 PetscInt i; 846 PetscErrorCode ierr; 847 848 PetscFunctionBegin; 849 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 850 PetscValidLogicalCollectiveInt(pc,n_is,2); 851 for (i=0;i<n_is;i++) { 852 PetscCheckSameComm(pc,1,ISForDofs[i],3); 853 PetscValidHeaderSpecific(ISForDofs[i],IS_CLASSID,3); 854 } 855 ierr = PetscTryMethod(pc,"PCBDDCSetDofsSplittingLocal_C",(PC,PetscInt,IS[]),(pc,n_is,ISForDofs));CHKERRQ(ierr); 856 PetscFunctionReturn(0); 857 } 858 /* -------------------------------------------------------------------------- */ 859 860 #undef __FUNCT__ 861 #define __FUNCT__ "PCBDDCSetDofsSplitting_BDDC" 862 static PetscErrorCode PCBDDCSetDofsSplitting_BDDC(PC pc,PetscInt n_is, IS ISForDofs[]) 863 { 864 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 865 PetscInt i; 866 PetscErrorCode ierr; 867 868 PetscFunctionBegin; 869 /* Destroy ISes if they were already set */ 870 for (i=0;i<pcbddc->n_ISForDofs;i++) { 871 ierr = ISDestroy(&pcbddc->ISForDofs[i]);CHKERRQ(ierr); 872 } 873 ierr = PetscFree(pcbddc->ISForDofs);CHKERRQ(ierr); 874 /* last user setting takes precendence -> destroy any other customization */ 875 for (i=0;i<pcbddc->n_ISForDofsLocal;i++) { 876 ierr = ISDestroy(&pcbddc->ISForDofsLocal[i]);CHKERRQ(ierr); 877 } 878 ierr = PetscFree(pcbddc->ISForDofsLocal);CHKERRQ(ierr); 879 pcbddc->n_ISForDofsLocal = 0; 880 /* allocate space then set */ 881 if (n_is) { 882 ierr = PetscMalloc1(n_is,&pcbddc->ISForDofs);CHKERRQ(ierr); 883 } 884 for (i=0;i<n_is;i++) { 885 ierr = PetscObjectReference((PetscObject)ISForDofs[i]);CHKERRQ(ierr); 886 pcbddc->ISForDofs[i]=ISForDofs[i]; 887 } 888 pcbddc->n_ISForDofs=n_is; 889 if (n_is) pcbddc->user_provided_isfordofs = PETSC_TRUE; 890 pcbddc->recompute_topography = PETSC_TRUE; 891 PetscFunctionReturn(0); 892 } 893 894 #undef __FUNCT__ 895 #define __FUNCT__ "PCBDDCSetDofsSplitting" 896 /*@ 897 PCBDDCSetDofsSplitting - Set index sets defining fields of the global matrix 898 899 Collective 900 901 Input Parameters: 902 + pc - the preconditioning context 903 . n_is - number of index sets defining the fields 904 - ISForDofs - array of IS describing the fields in global ordering 905 906 Level: intermediate 907 908 Notes: 909 Any process can list any global node. Not all nodes need to be listed: unlisted nodes will belong to the complement field. 910 911 .seealso: PCBDDC 912 @*/ 913 PetscErrorCode PCBDDCSetDofsSplitting(PC pc,PetscInt n_is, IS ISForDofs[]) 914 { 915 PetscInt i; 916 PetscErrorCode ierr; 917 918 PetscFunctionBegin; 919 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 920 PetscValidLogicalCollectiveInt(pc,n_is,2); 921 for (i=0;i<n_is;i++) { 922 PetscCheckSameComm(pc,1,ISForDofs[i],3); 923 PetscValidHeaderSpecific(ISForDofs[i],IS_CLASSID,3); 924 } 925 ierr = PetscTryMethod(pc,"PCBDDCSetDofsSplitting_C",(PC,PetscInt,IS[]),(pc,n_is,ISForDofs));CHKERRQ(ierr); 926 PetscFunctionReturn(0); 927 } 928 929 /* -------------------------------------------------------------------------- */ 930 #undef __FUNCT__ 931 #define __FUNCT__ "PCPreSolve_BDDC" 932 /* -------------------------------------------------------------------------- */ 933 /* 934 PCPreSolve_BDDC - Changes the right hand side and (if necessary) the initial 935 guess if a transformation of basis approach has been selected. 936 937 Input Parameter: 938 + pc - the preconditioner contex 939 940 Application Interface Routine: PCPreSolve() 941 942 Notes: 943 The interface routine PCPreSolve() is not usually called directly by 944 the user, but instead is called by KSPSolve(). 945 */ 946 static PetscErrorCode PCPreSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x) 947 { 948 PetscErrorCode ierr; 949 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 950 PC_IS *pcis = (PC_IS*)(pc->data); 951 Vec used_vec; 952 PetscBool copy_rhs = PETSC_TRUE; 953 954 PetscFunctionBegin; 955 /* if we are working with cg, one dirichlet solve can be avoided during Krylov iterations */ 956 if (ksp) { 957 PetscBool iscg; 958 ierr = PetscObjectTypeCompare((PetscObject)ksp,KSPCG,&iscg);CHKERRQ(ierr); 959 if (!iscg) { 960 ierr = PCBDDCSetUseExactDirichlet(pc,PETSC_FALSE);CHKERRQ(ierr); 961 } 962 } 963 /* Creates parallel work vectors used in presolve */ 964 if (!pcbddc->original_rhs) { 965 ierr = VecDuplicate(pcis->vec1_global,&pcbddc->original_rhs);CHKERRQ(ierr); 966 } 967 if (!pcbddc->temp_solution) { 968 ierr = VecDuplicate(pcis->vec1_global,&pcbddc->temp_solution);CHKERRQ(ierr); 969 } 970 971 if (x) { 972 ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr); 973 used_vec = x; 974 } else { /* it can only happen when calling PCBDDCMatFETIDPGetRHS */ 975 ierr = PetscObjectReference((PetscObject)pcbddc->temp_solution);CHKERRQ(ierr); 976 used_vec = pcbddc->temp_solution; 977 ierr = VecSet(used_vec,0.0);CHKERRQ(ierr); 978 } 979 980 /* hack into ksp data structure since PCPreSolve comes earlier than setting to zero the guess in src/ksp/ksp/interface/itfunc.c */ 981 if (ksp) { 982 /* store the flag for the initial guess since it will be restored back during PCPostSolve_BDDC */ 983 ierr = KSPGetInitialGuessNonzero(ksp,&pcbddc->ksp_guess_nonzero);CHKERRQ(ierr); 984 if (!pcbddc->ksp_guess_nonzero) { 985 ierr = VecSet(used_vec,0.0);CHKERRQ(ierr); 986 } 987 } 988 989 pcbddc->rhs_change = PETSC_FALSE; 990 /* Take into account zeroed rows -> change rhs and store solution removed */ 991 if (rhs) { 992 IS dirIS = NULL; 993 994 /* DirichletBoundariesLocal may not be consistent among neighbours; gets a dirichlet dofs IS from graph (may be cached) */ 995 ierr = PCBDDCGraphGetDirichletDofs(pcbddc->mat_graph,&dirIS);CHKERRQ(ierr); 996 if (dirIS) { 997 Mat_IS *matis = (Mat_IS*)pc->pmat->data; 998 PetscInt dirsize,i,*is_indices; 999 PetscScalar *array_x; 1000 const PetscScalar *array_diagonal; 1001 1002 ierr = MatGetDiagonal(pc->pmat,pcis->vec1_global);CHKERRQ(ierr); 1003 ierr = VecPointwiseDivide(pcis->vec1_global,rhs,pcis->vec1_global);CHKERRQ(ierr); 1004 ierr = VecScatterBegin(matis->rctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1005 ierr = VecScatterEnd(matis->rctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1006 ierr = VecScatterBegin(matis->rctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1007 ierr = VecScatterEnd(matis->rctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1008 ierr = ISGetLocalSize(dirIS,&dirsize);CHKERRQ(ierr); 1009 ierr = VecGetArray(pcis->vec1_N,&array_x);CHKERRQ(ierr); 1010 ierr = VecGetArrayRead(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr); 1011 ierr = ISGetIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr); 1012 for (i=0; i<dirsize; i++) array_x[is_indices[i]] = array_diagonal[is_indices[i]]; 1013 ierr = ISRestoreIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr); 1014 ierr = VecRestoreArrayRead(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr); 1015 ierr = VecRestoreArray(pcis->vec1_N,&array_x);CHKERRQ(ierr); 1016 ierr = VecScatterBegin(matis->rctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1017 ierr = VecScatterEnd(matis->rctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1018 pcbddc->rhs_change = PETSC_TRUE; 1019 ierr = ISDestroy(&dirIS);CHKERRQ(ierr); 1020 } 1021 } 1022 1023 /* remove the computed solution or the initial guess from the rhs */ 1024 if (pcbddc->rhs_change || (ksp && pcbddc->ksp_guess_nonzero) ) { 1025 /* store the original rhs */ 1026 if (copy_rhs) { 1027 ierr = VecCopy(rhs,pcbddc->original_rhs);CHKERRQ(ierr); 1028 copy_rhs = PETSC_FALSE; 1029 } 1030 pcbddc->rhs_change = PETSC_TRUE; 1031 ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr); 1032 ierr = MatMultAdd(pc->mat,used_vec,rhs,rhs);CHKERRQ(ierr); 1033 ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr); 1034 ierr = VecCopy(used_vec,pcbddc->temp_solution);CHKERRQ(ierr); 1035 if (ksp) { 1036 ierr = KSPSetInitialGuessNonzero(ksp,PETSC_FALSE);CHKERRQ(ierr); 1037 } 1038 } 1039 ierr = VecDestroy(&used_vec);CHKERRQ(ierr); 1040 1041 /* When using the benign trick: (TODO: what about FETI-DP?) 1042 - change rhs on pressures (iteration matrix is surely of type MATIS) 1043 - compute initial vector in benign space 1044 */ 1045 if (rhs && pcbddc->benign_saddle_point) { 1046 Mat_IS *matis = (Mat_IS*)(pc->mat->data); 1047 1048 /* store the original rhs */ 1049 if (copy_rhs) { 1050 ierr = VecCopy(rhs,pcbddc->original_rhs);CHKERRQ(ierr); 1051 copy_rhs = PETSC_FALSE; 1052 } 1053 1054 /* now change (locally) the basis */ 1055 ierr = VecScatterBegin(matis->rctx,rhs,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1056 ierr = VecScatterEnd(matis->rctx,rhs,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1057 if (pcbddc->benign_change) { 1058 ierr = MatMultTranspose(pcbddc->benign_change,pcis->vec1_N,pcis->vec2_N);CHKERRQ(ierr); 1059 ierr = VecScatterBegin(matis->rctx,pcis->vec2_N,rhs,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1060 ierr = VecScatterEnd(matis->rctx,pcis->vec2_N,rhs,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1061 /* change local iteration matrix */ 1062 ierr = MatPtAP(matis->A,pcbddc->benign_change,MAT_REUSE_MATRIX,2.0,&pcbddc->local_mat);CHKERRQ(ierr); 1063 ierr = MatDestroy(&pcbddc->benign_original_mat);CHKERRQ(ierr); 1064 ierr = PetscObjectReference((PetscObject)matis->A);CHKERRQ(ierr); 1065 pcbddc->benign_original_mat = matis->A; 1066 ierr = MatDestroy(&matis->A);CHKERRQ(ierr); 1067 ierr = PetscObjectReference((PetscObject)pcbddc->local_mat);CHKERRQ(ierr); 1068 matis->A = pcbddc->local_mat; 1069 ierr = MatDestroy(&pcbddc->local_mat);CHKERRQ(ierr); 1070 } else { 1071 ierr = VecScatterBegin(matis->rctx,pcis->vec1_N,rhs,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1072 ierr = VecScatterEnd(matis->rctx,pcis->vec1_N,rhs,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1073 } 1074 pcbddc->rhs_change = PETSC_TRUE; 1075 1076 /* compute u^*_h as in Xuemin Tu's thesis (see Section 4.8.1) */ 1077 /* TODO: what about Stokes? */ 1078 if (!pcbddc->benign_vec) { 1079 ierr = VecDuplicate(rhs,&pcbddc->benign_vec);CHKERRQ(ierr); 1080 } 1081 ierr = VecSet(pcis->vec1_global,0.);CHKERRQ(ierr); 1082 pcbddc->benign_p0 = 0.; 1083 if (pcbddc->benign_p0_lidx >= 0) { 1084 const PetscScalar *array; 1085 1086 ierr = VecGetArrayRead(pcis->vec2_N,&array);CHKERRQ(ierr); 1087 pcbddc->benign_p0 = array[pcbddc->benign_p0_lidx]; 1088 ierr = VecRestoreArrayRead(pcis->vec2_N,&array);CHKERRQ(ierr); 1089 } 1090 ierr = PCBDDCBenignGetOrSetP0(pc,pcis->vec1_global,PETSC_FALSE);CHKERRQ(ierr); 1091 ierr = PCApply_BDDC(pc,pcis->vec1_global,pcbddc->benign_vec);CHKERRQ(ierr); 1092 pcbddc->benign_p0 = 0.; 1093 ierr = PCBDDCBenignGetOrSetP0(pc,pcbddc->benign_vec,PETSC_FALSE);CHKERRQ(ierr); 1094 } 1095 1096 /* change rhs and iteration matrix if using the change of basis */ 1097 if (pcbddc->ChangeOfBasisMatrix) { 1098 PCBDDCChange_ctx change_ctx; 1099 1100 /* get change ctx */ 1101 ierr = MatShellGetContext(pcbddc->new_global_mat,&change_ctx);CHKERRQ(ierr); 1102 1103 /* set current iteration matrix inside change context (change of basis has been already set into the ctx during PCSetUp) */ 1104 ierr = MatDestroy(&change_ctx->original_mat);CHKERRQ(ierr); 1105 ierr = PetscObjectReference((PetscObject)pc->mat);CHKERRQ(ierr); 1106 change_ctx->original_mat = pc->mat; 1107 1108 /* change iteration matrix */ 1109 ierr = MatDestroy(&pc->mat);CHKERRQ(ierr); 1110 ierr = PetscObjectReference((PetscObject)pcbddc->new_global_mat);CHKERRQ(ierr); 1111 pc->mat = pcbddc->new_global_mat; 1112 1113 /* store the original rhs */ 1114 if (copy_rhs) { 1115 ierr = VecCopy(rhs,pcbddc->original_rhs);CHKERRQ(ierr); 1116 copy_rhs = PETSC_FALSE; 1117 } 1118 1119 /* change rhs */ 1120 ierr = MatMultTranspose(change_ctx->global_change,rhs,pcis->vec1_global);CHKERRQ(ierr); 1121 ierr = VecCopy(pcis->vec1_global,rhs);CHKERRQ(ierr); 1122 pcbddc->rhs_change = PETSC_TRUE; 1123 } 1124 1125 /* remove non-benign solution from the rhs */ 1126 if (pcbddc->benign_saddle_point) { 1127 /* store the original rhs */ 1128 if (copy_rhs) { 1129 ierr = VecCopy(rhs,pcbddc->original_rhs);CHKERRQ(ierr); 1130 copy_rhs = PETSC_FALSE; 1131 } 1132 ierr = VecScale(pcbddc->benign_vec,-1.0);CHKERRQ(ierr); 1133 ierr = MatMultAdd(pc->mat,pcbddc->benign_vec,rhs,rhs);CHKERRQ(ierr); 1134 pcbddc->rhs_change = PETSC_TRUE; 1135 } 1136 1137 /* set initial guess if using PCG */ 1138 if (x && pcbddc->use_exact_dirichlet_trick) { 1139 ierr = VecSet(x,0.0);CHKERRQ(ierr); 1140 ierr = VecScatterBegin(pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1141 ierr = VecScatterEnd(pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1142 ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr); 1143 ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1144 ierr = VecScatterEnd(pcis->global_to_D,pcis->vec2_D,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1145 if (ksp) { 1146 ierr = KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);CHKERRQ(ierr); 1147 } 1148 } 1149 1150 /* remove nullspace if present */ 1151 if (ksp && x && pcbddc->NullSpace) { 1152 ierr = MatNullSpaceRemove(pcbddc->NullSpace,x);CHKERRQ(ierr); 1153 /* store the original rhs */ 1154 if (copy_rhs) { 1155 ierr = VecCopy(rhs,pcbddc->original_rhs);CHKERRQ(ierr); 1156 copy_rhs = PETSC_FALSE; 1157 } 1158 pcbddc->rhs_change = PETSC_TRUE; 1159 ierr = MatNullSpaceRemove(pcbddc->NullSpace,rhs);CHKERRQ(ierr); 1160 } 1161 PetscFunctionReturn(0); 1162 } 1163 1164 /* -------------------------------------------------------------------------- */ 1165 #undef __FUNCT__ 1166 #define __FUNCT__ "PCPostSolve_BDDC" 1167 /* -------------------------------------------------------------------------- */ 1168 /* 1169 PCPostSolve_BDDC - Changes the computed solution if a transformation of basis 1170 approach has been selected. Also, restores rhs to its original state. 1171 1172 Input Parameter: 1173 + pc - the preconditioner contex 1174 1175 Application Interface Routine: PCPostSolve() 1176 1177 Notes: 1178 The interface routine PCPostSolve() is not usually called directly by 1179 the user, but instead is called by KSPSolve(). 1180 */ 1181 static PetscErrorCode PCPostSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x) 1182 { 1183 PetscErrorCode ierr; 1184 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 1185 1186 PetscFunctionBegin; 1187 if (pcbddc->ChangeOfBasisMatrix) { 1188 PCBDDCChange_ctx change_ctx; 1189 1190 /* get change ctx */ 1191 ierr = MatShellGetContext(pcbddc->new_global_mat,&change_ctx);CHKERRQ(ierr); 1192 1193 /* restore iteration matrix */ 1194 ierr = MatDestroy(&pc->mat);CHKERRQ(ierr); 1195 ierr = PetscObjectReference((PetscObject)change_ctx->original_mat);CHKERRQ(ierr); 1196 pc->mat = change_ctx->original_mat; 1197 1198 /* need to restore the local matrices */ 1199 if (pcbddc->benign_saddle_point && pcbddc->benign_change) { 1200 Mat_IS *matis = (Mat_IS*)pc->mat->data; 1201 1202 ierr = MatDestroy(&matis->A);CHKERRQ(ierr); 1203 ierr = PetscObjectReference((PetscObject)pcbddc->benign_original_mat);CHKERRQ(ierr); 1204 matis->A = pcbddc->benign_original_mat; 1205 ierr = MatDestroy(&pcbddc->benign_original_mat);CHKERRQ(ierr); 1206 } 1207 1208 /* get solution in original basis */ 1209 if (x) { 1210 PC_IS *pcis = (PC_IS*)(pc->data); 1211 1212 /* restore solution on pressures */ 1213 if (pcbddc->benign_saddle_point) { 1214 Mat_IS *matis = (Mat_IS*)pc->mat->data; 1215 1216 /* add non-benign solution */ 1217 ierr = VecAXPY(x,-1.0,pcbddc->benign_vec);CHKERRQ(ierr); 1218 1219 /* change basis on pressures for x */ 1220 ierr = VecScatterBegin(matis->rctx,x,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1221 ierr = VecScatterEnd(matis->rctx,x,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1222 if (pcbddc->benign_change) { 1223 1224 ierr = MatMult(pcbddc->benign_change,pcis->vec1_N,pcis->vec2_N);CHKERRQ(ierr); 1225 ierr = VecScatterBegin(matis->rctx,pcis->vec2_N,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1226 ierr = VecScatterEnd(matis->rctx,pcis->vec2_N,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1227 } else { 1228 ierr = VecScatterBegin(matis->rctx,pcis->vec1_N,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1229 ierr = VecScatterEnd(matis->rctx,pcis->vec1_N,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1230 } 1231 } 1232 /* change basis on x */ 1233 ierr = MatMult(change_ctx->global_change,x,pcis->vec1_global);CHKERRQ(ierr); 1234 ierr = VecCopy(pcis->vec1_global,x);CHKERRQ(ierr); 1235 } 1236 } 1237 1238 /* add solution removed in presolve */ 1239 if (x && pcbddc->rhs_change) { 1240 ierr = VecAXPY(x,1.0,pcbddc->temp_solution);CHKERRQ(ierr); 1241 } 1242 1243 /* restore rhs to its original state */ 1244 if (rhs && pcbddc->rhs_change) { 1245 ierr = VecCopy(pcbddc->original_rhs,rhs);CHKERRQ(ierr); 1246 } 1247 pcbddc->rhs_change = PETSC_FALSE; 1248 1249 /* restore ksp guess state */ 1250 if (ksp) { 1251 ierr = KSPSetInitialGuessNonzero(ksp,pcbddc->ksp_guess_nonzero);CHKERRQ(ierr); 1252 } 1253 PetscFunctionReturn(0); 1254 } 1255 /* -------------------------------------------------------------------------- */ 1256 #undef __FUNCT__ 1257 #define __FUNCT__ "PCSetUp_BDDC" 1258 /* -------------------------------------------------------------------------- */ 1259 /* 1260 PCSetUp_BDDC - Prepares for the use of the BDDC preconditioner 1261 by setting data structures and options. 1262 1263 Input Parameter: 1264 + pc - the preconditioner context 1265 1266 Application Interface Routine: PCSetUp() 1267 1268 Notes: 1269 The interface routine PCSetUp() is not usually called directly by 1270 the user, but instead is called by PCApply() if necessary. 1271 */ 1272 PetscErrorCode PCSetUp_BDDC(PC pc) 1273 { 1274 PetscErrorCode ierr; 1275 PC_BDDC* pcbddc = (PC_BDDC*)pc->data; 1276 Mat_IS* matis; 1277 MatNullSpace nearnullspace; 1278 IS zerodiag = NULL; 1279 PetscInt nrows,ncols; 1280 PetscBool computetopography,computesolvers,computesubschurs; 1281 PetscBool computeconstraintsmatrix; 1282 PetscBool new_nearnullspace_provided,ismatis; 1283 1284 PetscFunctionBegin; 1285 ierr = PetscObjectTypeCompare((PetscObject)pc->pmat,MATIS,&ismatis);CHKERRQ(ierr); 1286 if (!ismatis) { 1287 SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"PCBDDC preconditioner requires matrix of type MATIS"); 1288 } 1289 ierr = MatGetSize(pc->pmat,&nrows,&ncols);CHKERRQ(ierr); 1290 if (nrows != ncols) { 1291 SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"PCBDDC preconditioner requires a square preconditioning matrix"); 1292 } 1293 matis = (Mat_IS*)pc->pmat->data; 1294 /* the following lines of code should be replaced by a better logic between PCIS, PCNN, PCBDDC and other future nonoverlapping preconditioners */ 1295 /* For BDDC we need to define a local "Neumann" problem different to that defined in PCISSetup 1296 Also, BDDC directly build the Dirichlet problem */ 1297 /* split work */ 1298 if (pc->setupcalled) { 1299 if (pc->flag == SAME_NONZERO_PATTERN) { 1300 computetopography = PETSC_FALSE; 1301 computesolvers = PETSC_TRUE; 1302 } else { /* DIFFERENT_NONZERO_PATTERN */ 1303 computetopography = PETSC_TRUE; 1304 computesolvers = PETSC_TRUE; 1305 } 1306 } else { 1307 computetopography = PETSC_TRUE; 1308 computesolvers = PETSC_TRUE; 1309 } 1310 if (pcbddc->recompute_topography) { 1311 computetopography = PETSC_TRUE; 1312 } 1313 computeconstraintsmatrix = PETSC_FALSE; 1314 1315 /* check parameters' compatibility */ 1316 if (pcbddc->adaptive_threshold > 0.0 && !pcbddc->use_deluxe_scaling) { 1317 SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Cannot compute adaptive constraints without deluxe scaling. Rerun with -pc_bddc_use_deluxe_scaling"); 1318 } 1319 pcbddc->adaptive_selection = (PetscBool)(pcbddc->adaptive_threshold > 0.0 && pcbddc->use_deluxe_scaling); 1320 if (pcbddc->adaptive_selection) pcbddc->use_faces = PETSC_TRUE; 1321 1322 computesubschurs = (PetscBool)(pcbddc->adaptive_selection || pcbddc->use_deluxe_scaling); 1323 if (pcbddc->faster_deluxe && pcbddc->adaptive_selection && pcbddc->use_change_of_basis) { 1324 SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Cannot compute faster deluxe if adaptivity and change of basis are both requested. Rerun with -pc_bddc_deluxe_faster false"); 1325 } 1326 1327 /* check if the iteration matrix is of type MATIS in case the benign trick has been requested */ 1328 ierr = PetscObjectTypeCompare((PetscObject)pc->mat,MATIS,&ismatis);CHKERRQ(ierr); 1329 if (pcbddc->benign_saddle_point && !ismatis) { 1330 SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"PCBDDC preconditioner with benign subspace trick requires the iteration matrix to be of type MATIS"); 1331 } 1332 1333 /* Get stdout for dbg */ 1334 if (pcbddc->dbg_flag) { 1335 if (!pcbddc->dbg_viewer) { 1336 pcbddc->dbg_viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)pc)); 1337 ierr = PetscViewerASCIISynchronizedAllow(pcbddc->dbg_viewer,PETSC_TRUE);CHKERRQ(ierr); 1338 } 1339 ierr = PetscViewerASCIIAddTab(pcbddc->dbg_viewer,2*pcbddc->current_level);CHKERRQ(ierr); 1340 } 1341 1342 if (pcbddc->user_ChangeOfBasisMatrix) { 1343 /* use_change_of_basis flag is used to automatically compute a change of basis from constraints */ 1344 pcbddc->use_change_of_basis = PETSC_FALSE; 1345 ierr = PCBDDCComputeLocalMatrix(pc,pcbddc->user_ChangeOfBasisMatrix);CHKERRQ(ierr); 1346 } else { 1347 ierr = MatDestroy(&pcbddc->local_mat);CHKERRQ(ierr); 1348 ierr = PetscObjectReference((PetscObject)matis->A);CHKERRQ(ierr); 1349 pcbddc->local_mat = matis->A; 1350 } 1351 1352 /* change basis on local pressures (aka zerodiag dofs) */ 1353 if (pcbddc->benign_saddle_point) { 1354 /* detect local saddle point and change the basis in pcbddc->local_mat */ 1355 ierr = PCBDDCBenignDetectSaddlePoint(pc,&zerodiag);CHKERRQ(ierr); 1356 /* pop B0 mat from pcbddc->local_mat */ 1357 ierr = PCBDDCBenignPopOrPushB0(pc,PETSC_TRUE);CHKERRQ(ierr); 1358 } 1359 1360 /* propagate relevant information */ 1361 /* workaround for reals */ 1362 #if !defined(PETSC_USE_COMPLEX) 1363 if (matis->A->symmetric_set) { 1364 ierr = MatSetOption(pcbddc->local_mat,MAT_HERMITIAN,matis->A->symmetric);CHKERRQ(ierr); 1365 } 1366 #endif 1367 if (matis->A->symmetric_set) { 1368 ierr = MatSetOption(pcbddc->local_mat,MAT_SYMMETRIC,matis->A->symmetric);CHKERRQ(ierr); 1369 } 1370 if (matis->A->spd_set) { 1371 ierr = MatSetOption(pcbddc->local_mat,MAT_SPD,matis->A->spd);CHKERRQ(ierr); 1372 } 1373 1374 /* Set up all the "iterative substructuring" common block without computing solvers */ 1375 { 1376 Mat temp_mat; 1377 1378 temp_mat = matis->A; 1379 matis->A = pcbddc->local_mat; 1380 ierr = PCISSetUp(pc,PETSC_FALSE);CHKERRQ(ierr); 1381 pcbddc->local_mat = matis->A; 1382 matis->A = temp_mat; 1383 } 1384 1385 /* Analyze interface */ 1386 if (computetopography) { 1387 ierr = PCBDDCAnalyzeInterface(pc);CHKERRQ(ierr); 1388 computeconstraintsmatrix = PETSC_TRUE; 1389 } 1390 1391 /* check existence of a divergence free extension, i.e. 1392 b(v_I,p_0) = 0 for all v_I (raise error if not). 1393 Also, check that PCBDDCBenignGetOrSetP0 works */ 1394 #if defined(PETSC_USE_DEBUG) 1395 if (pcbddc->benign_saddle_point) { 1396 ierr = PCBDDCBenignCheck(pc,zerodiag);CHKERRQ(ierr); 1397 } 1398 #endif 1399 1400 /* Setup local dirichlet solver ksp_D and sub_schurs solvers */ 1401 if (computesolvers) { 1402 PCBDDCSubSchurs sub_schurs=pcbddc->sub_schurs; 1403 1404 if (computesubschurs && computetopography) { 1405 ierr = PCBDDCInitSubSchurs(pc);CHKERRQ(ierr); 1406 } 1407 if (sub_schurs->use_mumps) { 1408 if (computesubschurs) { 1409 ierr = PCBDDCSetUpSubSchurs(pc);CHKERRQ(ierr); 1410 } 1411 ierr = PCBDDCSetUpLocalSolvers(pc,PETSC_TRUE,PETSC_FALSE);CHKERRQ(ierr); 1412 } else { 1413 ierr = PCBDDCSetUpLocalSolvers(pc,PETSC_TRUE,PETSC_FALSE);CHKERRQ(ierr); 1414 if (computesubschurs) { 1415 ierr = PCBDDCSetUpSubSchurs(pc);CHKERRQ(ierr); 1416 } 1417 } 1418 if (pcbddc->adaptive_selection) { 1419 ierr = PCBDDCAdaptiveSelection(pc);CHKERRQ(ierr); 1420 computeconstraintsmatrix = PETSC_TRUE; 1421 } 1422 } 1423 1424 /* infer if NullSpace object attached to Mat via MatSetNearNullSpace has changed */ 1425 new_nearnullspace_provided = PETSC_FALSE; 1426 ierr = MatGetNearNullSpace(pc->pmat,&nearnullspace);CHKERRQ(ierr); 1427 if (pcbddc->onearnullspace) { /* already used nearnullspace */ 1428 if (!nearnullspace) { /* near null space attached to mat has been destroyed */ 1429 new_nearnullspace_provided = PETSC_TRUE; 1430 } else { 1431 /* determine if the two nullspaces are different (should be lightweight) */ 1432 if (nearnullspace != pcbddc->onearnullspace) { 1433 new_nearnullspace_provided = PETSC_TRUE; 1434 } else { /* maybe the user has changed the content of the nearnullspace so check vectors ObjectStateId */ 1435 PetscInt i; 1436 const Vec *nearnullvecs; 1437 PetscObjectState state; 1438 PetscInt nnsp_size; 1439 ierr = MatNullSpaceGetVecs(nearnullspace,NULL,&nnsp_size,&nearnullvecs);CHKERRQ(ierr); 1440 for (i=0;i<nnsp_size;i++) { 1441 ierr = PetscObjectStateGet((PetscObject)nearnullvecs[i],&state);CHKERRQ(ierr); 1442 if (pcbddc->onearnullvecs_state[i] != state) { 1443 new_nearnullspace_provided = PETSC_TRUE; 1444 break; 1445 } 1446 } 1447 } 1448 } 1449 } else { 1450 if (!nearnullspace) { /* both nearnullspaces are null */ 1451 new_nearnullspace_provided = PETSC_FALSE; 1452 } else { /* nearnullspace attached later */ 1453 new_nearnullspace_provided = PETSC_TRUE; 1454 } 1455 } 1456 1457 /* Setup constraints and related work vectors */ 1458 /* reset primal space flags */ 1459 pcbddc->new_primal_space = PETSC_FALSE; 1460 pcbddc->new_primal_space_local = PETSC_FALSE; 1461 if (computeconstraintsmatrix || new_nearnullspace_provided) { 1462 /* It also sets the primal space flags */ 1463 ierr = PCBDDCConstraintsSetUp(pc);CHKERRQ(ierr); 1464 /* Allocate needed local vectors (which depends on quantities defined during ConstraintsSetUp) */ 1465 ierr = PCBDDCSetUpLocalWorkVectors(pc);CHKERRQ(ierr); 1466 } 1467 1468 if (computesolvers || pcbddc->new_primal_space) { 1469 if (pcbddc->use_change_of_basis) { 1470 PC_IS *pcis = (PC_IS*)(pc->data); 1471 Mat temp_mat = NULL; 1472 1473 if (zerodiag) { 1474 /* insert B0 in pcbddc->local_mat */ 1475 ierr = PCBDDCBenignPopOrPushB0(pc,PETSC_FALSE);CHKERRQ(ierr); 1476 /* swap local matrices */ 1477 ierr = MatISGetLocalMat(pc->pmat,&temp_mat);CHKERRQ(ierr); 1478 ierr = PetscObjectReference((PetscObject)temp_mat);CHKERRQ(ierr); 1479 ierr = MatISSetLocalMat(pc->pmat,pcbddc->local_mat);CHKERRQ(ierr); 1480 ierr = MatDestroy(&pcbddc->local_mat);CHKERRQ(ierr); 1481 } 1482 ierr = PCBDDCComputeLocalMatrix(pc,pcbddc->ChangeOfBasisMatrix);CHKERRQ(ierr); 1483 if (zerodiag) { 1484 /* restore original matrix */ 1485 ierr = MatISSetLocalMat(pc->pmat,temp_mat);CHKERRQ(ierr); 1486 ierr = PetscObjectDereference((PetscObject)temp_mat);CHKERRQ(ierr); 1487 /* pop B0 from pcbddc->local_mat */ 1488 ierr = PCBDDCBenignPopOrPushB0(pc,PETSC_TRUE);CHKERRQ(ierr); 1489 } 1490 /* get submatrices */ 1491 ierr = MatDestroy(&pcis->A_IB);CHKERRQ(ierr); 1492 ierr = MatDestroy(&pcis->A_BI);CHKERRQ(ierr); 1493 ierr = MatDestroy(&pcis->A_BB);CHKERRQ(ierr); 1494 ierr = MatGetSubMatrix(pcbddc->local_mat,pcis->is_B_local,pcis->is_B_local,MAT_INITIAL_MATRIX,&pcis->A_BB);CHKERRQ(ierr); 1495 ierr = MatGetSubMatrix(pcbddc->local_mat,pcis->is_I_local,pcis->is_B_local,MAT_INITIAL_MATRIX,&pcis->A_IB);CHKERRQ(ierr); 1496 ierr = MatGetSubMatrix(pcbddc->local_mat,pcis->is_B_local,pcis->is_I_local,MAT_INITIAL_MATRIX,&pcis->A_BI);CHKERRQ(ierr); 1497 /* set flag in pcis to not reuse submatrices during PCISCreate */ 1498 pcis->reusesubmatrices = PETSC_FALSE; 1499 } else if (!pcbddc->user_ChangeOfBasisMatrix) { 1500 ierr = MatDestroy(&pcbddc->local_mat);CHKERRQ(ierr); 1501 ierr = PetscObjectReference((PetscObject)matis->A);CHKERRQ(ierr); 1502 pcbddc->local_mat = matis->A; 1503 } 1504 /* SetUp coarse and local Neumann solvers */ 1505 ierr = PCBDDCSetUpSolvers(pc);CHKERRQ(ierr); 1506 /* SetUp Scaling operator */ 1507 ierr = PCBDDCScalingSetUp(pc);CHKERRQ(ierr); 1508 } 1509 ierr = ISDestroy(&zerodiag);CHKERRQ(ierr); 1510 1511 if (pcbddc->dbg_flag) { 1512 ierr = PetscViewerASCIISubtractTab(pcbddc->dbg_viewer,2*pcbddc->current_level);CHKERRQ(ierr); 1513 } 1514 PetscFunctionReturn(0); 1515 } 1516 1517 /* -------------------------------------------------------------------------- */ 1518 /* 1519 PCApply_BDDC - Applies the BDDC operator to a vector. 1520 1521 Input Parameters: 1522 + pc - the preconditioner context 1523 - r - input vector (global) 1524 1525 Output Parameter: 1526 . z - output vector (global) 1527 1528 Application Interface Routine: PCApply() 1529 */ 1530 #undef __FUNCT__ 1531 #define __FUNCT__ "PCApply_BDDC" 1532 PetscErrorCode PCApply_BDDC(PC pc,Vec r,Vec z) 1533 { 1534 PC_IS *pcis = (PC_IS*)(pc->data); 1535 PC_BDDC *pcbddc = (PC_BDDC*)(pc->data); 1536 PetscInt n_B = pcis->n_B, n_D = pcis->n - n_B; 1537 PetscErrorCode ierr; 1538 const PetscScalar one = 1.0; 1539 const PetscScalar m_one = -1.0; 1540 const PetscScalar zero = 0.0; 1541 1542 /* This code is similar to that provided in nn.c for PCNN 1543 NN interface preconditioner changed to BDDC 1544 Added support for M_3 preconditioner in the reference article (code is active if pcbddc->switch_static == PETSC_TRUE) */ 1545 1546 PetscFunctionBegin; 1547 if (pcbddc->benign_saddle_point) { /* get p0 from r */ 1548 ierr = PCBDDCBenignGetOrSetP0(pc,r,PETSC_TRUE);CHKERRQ(ierr); 1549 } 1550 if (!pcbddc->use_exact_dirichlet_trick) { 1551 ierr = VecCopy(r,z);CHKERRQ(ierr); 1552 /* First Dirichlet solve */ 1553 ierr = VecScatterBegin(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1554 ierr = VecScatterEnd(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1555 /* 1556 Assembling right hand side for BDDC operator 1557 - pcis->vec1_D for the Dirichlet part (if needed, i.e. pcbddc->switch_static == PETSC_TRUE) 1558 - pcis->vec1_B the interface part of the global vector z 1559 */ 1560 if (n_D) { 1561 ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr); 1562 ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr); 1563 if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec2_D,pcis->vec1_D,pcis->vec1_D);CHKERRQ(ierr); } 1564 ierr = MatMult(pcis->A_BI,pcis->vec2_D,pcis->vec1_B);CHKERRQ(ierr); 1565 } else { 1566 ierr = VecSet(pcis->vec1_B,zero);CHKERRQ(ierr); 1567 } 1568 ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1569 ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1570 ierr = PCBDDCScalingRestriction(pc,z,pcis->vec1_B);CHKERRQ(ierr); 1571 } else { 1572 if (pcbddc->switch_static) { 1573 ierr = VecSet(pcis->vec1_D,zero);CHKERRQ(ierr); 1574 } 1575 ierr = PCBDDCScalingRestriction(pc,r,pcis->vec1_B);CHKERRQ(ierr); 1576 } 1577 1578 /* Apply interface preconditioner 1579 input/output vecs: pcis->vec1_B and pcis->vec1_D */ 1580 ierr = PCBDDCApplyInterfacePreconditioner(pc,PETSC_FALSE);CHKERRQ(ierr); 1581 1582 /* Apply transpose of partition of unity operator */ 1583 ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,z);CHKERRQ(ierr); 1584 1585 /* Second Dirichlet solve and assembling of output */ 1586 ierr = VecScatterBegin(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1587 ierr = VecScatterEnd(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1588 if (n_B) { 1589 ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec3_D);CHKERRQ(ierr); 1590 if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec1_D,pcis->vec3_D,pcis->vec3_D);CHKERRQ(ierr); } 1591 } else if (pcbddc->switch_static) { 1592 ierr = MatMult(pcis->A_II,pcis->vec1_D,pcis->vec3_D);CHKERRQ(ierr); 1593 } 1594 ierr = KSPSolve(pcbddc->ksp_D,pcis->vec3_D,pcis->vec4_D);CHKERRQ(ierr); 1595 1596 if (!pcbddc->use_exact_dirichlet_trick) { 1597 if (pcbddc->switch_static) { 1598 ierr = VecAXPBYPCZ(pcis->vec2_D,m_one,one,m_one,pcis->vec4_D,pcis->vec1_D);CHKERRQ(ierr); 1599 } else { 1600 ierr = VecAXPBY(pcis->vec2_D,m_one,m_one,pcis->vec4_D);CHKERRQ(ierr); 1601 } 1602 ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1603 ierr = VecScatterEnd(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1604 } else { 1605 if (pcbddc->switch_static) { 1606 ierr = VecAXPBY(pcis->vec4_D,one,m_one,pcis->vec1_D);CHKERRQ(ierr); 1607 } else { 1608 ierr = VecScale(pcis->vec4_D,m_one);CHKERRQ(ierr); 1609 } 1610 ierr = VecScatterBegin(pcis->global_to_D,pcis->vec4_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1611 ierr = VecScatterEnd(pcis->global_to_D,pcis->vec4_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1612 } 1613 1614 if (pcbddc->benign_saddle_point) { /* set p0 (computed in PCBDDCApplyInterface) */ 1615 ierr = PCBDDCBenignGetOrSetP0(pc,z,PETSC_FALSE);CHKERRQ(ierr); 1616 } 1617 PetscFunctionReturn(0); 1618 } 1619 1620 /* -------------------------------------------------------------------------- */ 1621 /* 1622 PCApplyTranspose_BDDC - Applies the transpose of the BDDC operator to a vector. 1623 1624 Input Parameters: 1625 + pc - the preconditioner context 1626 - r - input vector (global) 1627 1628 Output Parameter: 1629 . z - output vector (global) 1630 1631 Application Interface Routine: PCApplyTranspose() 1632 */ 1633 #undef __FUNCT__ 1634 #define __FUNCT__ "PCApplyTranspose_BDDC" 1635 PetscErrorCode PCApplyTranspose_BDDC(PC pc,Vec r,Vec z) 1636 { 1637 PC_IS *pcis = (PC_IS*)(pc->data); 1638 PC_BDDC *pcbddc = (PC_BDDC*)(pc->data); 1639 PetscInt n_B = pcis->n_B, n_D = pcis->n - n_B; 1640 PetscErrorCode ierr; 1641 const PetscScalar one = 1.0; 1642 const PetscScalar m_one = -1.0; 1643 const PetscScalar zero = 0.0; 1644 1645 PetscFunctionBegin; 1646 if (pcbddc->benign_saddle_point) { /* get p0 from r */ 1647 ierr = PCBDDCBenignGetOrSetP0(pc,r,PETSC_TRUE);CHKERRQ(ierr); 1648 } 1649 if (!pcbddc->use_exact_dirichlet_trick) { 1650 ierr = VecCopy(r,z);CHKERRQ(ierr); 1651 /* First Dirichlet solve */ 1652 ierr = VecScatterBegin(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1653 ierr = VecScatterEnd(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1654 /* 1655 Assembling right hand side for BDDC operator 1656 - pcis->vec1_D for the Dirichlet part (if needed, i.e. pcbddc->switch_static == PETSC_TRUE) 1657 - pcis->vec1_B the interface part of the global vector z 1658 */ 1659 if (n_D) { 1660 ierr = KSPSolveTranspose(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr); 1661 ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr); 1662 if (pcbddc->switch_static) { ierr = MatMultTransposeAdd(pcis->A_II,pcis->vec2_D,pcis->vec1_D,pcis->vec1_D);CHKERRQ(ierr); } 1663 ierr = MatMultTranspose(pcis->A_IB,pcis->vec2_D,pcis->vec1_B);CHKERRQ(ierr); 1664 } else { 1665 ierr = VecSet(pcis->vec1_B,zero);CHKERRQ(ierr); 1666 } 1667 ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1668 ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1669 ierr = PCBDDCScalingRestriction(pc,z,pcis->vec1_B);CHKERRQ(ierr); 1670 } else { 1671 if (pcbddc->switch_static) { 1672 ierr = VecSet(pcis->vec1_D,zero);CHKERRQ(ierr); 1673 } 1674 ierr = PCBDDCScalingRestriction(pc,r,pcis->vec1_B);CHKERRQ(ierr); 1675 } 1676 1677 /* Apply interface preconditioner 1678 input/output vecs: pcis->vec1_B and pcis->vec1_D */ 1679 ierr = PCBDDCApplyInterfacePreconditioner(pc,PETSC_TRUE);CHKERRQ(ierr); 1680 1681 /* Apply transpose of partition of unity operator */ 1682 ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,z);CHKERRQ(ierr); 1683 1684 /* Second Dirichlet solve and assembling of output */ 1685 ierr = VecScatterBegin(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1686 ierr = VecScatterEnd(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1687 if (n_B) { 1688 ierr = MatMultTranspose(pcis->A_BI,pcis->vec1_B,pcis->vec3_D);CHKERRQ(ierr); 1689 if (pcbddc->switch_static) { ierr = MatMultTransposeAdd(pcis->A_II,pcis->vec1_D,pcis->vec3_D,pcis->vec3_D);CHKERRQ(ierr); } 1690 } else if (pcbddc->switch_static) { 1691 ierr = MatMultTranspose(pcis->A_II,pcis->vec1_D,pcis->vec3_D);CHKERRQ(ierr); 1692 } 1693 ierr = KSPSolveTranspose(pcbddc->ksp_D,pcis->vec3_D,pcis->vec4_D);CHKERRQ(ierr); 1694 if (!pcbddc->use_exact_dirichlet_trick) { 1695 if (pcbddc->switch_static) { 1696 ierr = VecAXPBYPCZ(pcis->vec2_D,m_one,one,m_one,pcis->vec4_D,pcis->vec1_D);CHKERRQ(ierr); 1697 } else { 1698 ierr = VecAXPBY(pcis->vec2_D,m_one,m_one,pcis->vec4_D);CHKERRQ(ierr); 1699 } 1700 ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1701 ierr = VecScatterEnd(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1702 } else { 1703 if (pcbddc->switch_static) { 1704 ierr = VecAXPBY(pcis->vec4_D,one,m_one,pcis->vec1_D);CHKERRQ(ierr); 1705 } else { 1706 ierr = VecScale(pcis->vec4_D,m_one);CHKERRQ(ierr); 1707 } 1708 ierr = VecScatterBegin(pcis->global_to_D,pcis->vec4_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1709 ierr = VecScatterEnd(pcis->global_to_D,pcis->vec4_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1710 } 1711 if (pcbddc->benign_saddle_point) { /* set p0 (computed in PCBDDCApplyInterface) */ 1712 ierr = PCBDDCBenignGetOrSetP0(pc,z,PETSC_FALSE);CHKERRQ(ierr); 1713 } 1714 PetscFunctionReturn(0); 1715 } 1716 /* -------------------------------------------------------------------------- */ 1717 1718 #undef __FUNCT__ 1719 #define __FUNCT__ "PCDestroy_BDDC" 1720 PetscErrorCode PCDestroy_BDDC(PC pc) 1721 { 1722 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 1723 PetscErrorCode ierr; 1724 1725 PetscFunctionBegin; 1726 /* free data created by PCIS */ 1727 ierr = PCISDestroy(pc);CHKERRQ(ierr); 1728 /* free BDDC custom data */ 1729 ierr = PCBDDCResetCustomization(pc);CHKERRQ(ierr); 1730 /* destroy objects related to topography */ 1731 ierr = PCBDDCResetTopography(pc);CHKERRQ(ierr); 1732 /* free allocated graph structure */ 1733 ierr = PetscFree(pcbddc->mat_graph);CHKERRQ(ierr); 1734 /* free allocated sub schurs structure */ 1735 ierr = PetscFree(pcbddc->sub_schurs);CHKERRQ(ierr); 1736 /* destroy objects for scaling operator */ 1737 ierr = PCBDDCScalingDestroy(pc);CHKERRQ(ierr); 1738 ierr = PetscFree(pcbddc->deluxe_ctx);CHKERRQ(ierr); 1739 /* free solvers stuff */ 1740 ierr = PCBDDCResetSolvers(pc);CHKERRQ(ierr); 1741 /* free global vectors needed in presolve */ 1742 ierr = VecDestroy(&pcbddc->temp_solution);CHKERRQ(ierr); 1743 ierr = VecDestroy(&pcbddc->original_rhs);CHKERRQ(ierr); 1744 /* free stuff for change of basis hooks */ 1745 if (pcbddc->new_global_mat) { 1746 PCBDDCChange_ctx change_ctx; 1747 ierr = MatShellGetContext(pcbddc->new_global_mat,&change_ctx);CHKERRQ(ierr); 1748 ierr = MatDestroy(&change_ctx->original_mat);CHKERRQ(ierr); 1749 ierr = MatDestroy(&change_ctx->global_change);CHKERRQ(ierr); 1750 ierr = VecDestroyVecs(2,&change_ctx->work);CHKERRQ(ierr); 1751 ierr = PetscFree(change_ctx);CHKERRQ(ierr); 1752 } 1753 ierr = MatDestroy(&pcbddc->new_global_mat);CHKERRQ(ierr); 1754 /* remove functions */ 1755 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetChangeOfBasisMat_C",NULL);CHKERRQ(ierr); 1756 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",NULL);CHKERRQ(ierr); 1757 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",NULL);CHKERRQ(ierr); 1758 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",NULL);CHKERRQ(ierr); 1759 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",NULL);CHKERRQ(ierr); 1760 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",NULL);CHKERRQ(ierr); 1761 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",NULL);CHKERRQ(ierr); 1762 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",NULL);CHKERRQ(ierr); 1763 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundariesLocal_C",NULL);CHKERRQ(ierr); 1764 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",NULL);CHKERRQ(ierr); 1765 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundariesLocal_C",NULL);CHKERRQ(ierr); 1766 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",NULL);CHKERRQ(ierr); 1767 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundariesLocal_C",NULL);CHKERRQ(ierr); 1768 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",NULL);CHKERRQ(ierr); 1769 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundariesLocal_C",NULL);CHKERRQ(ierr); 1770 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",NULL);CHKERRQ(ierr); 1771 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplittingLocal_C",NULL);CHKERRQ(ierr); 1772 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",NULL);CHKERRQ(ierr); 1773 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",NULL);CHKERRQ(ierr); 1774 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",NULL);CHKERRQ(ierr); 1775 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",NULL);CHKERRQ(ierr); 1776 /* Free the private data structure */ 1777 ierr = PetscFree(pc->data);CHKERRQ(ierr); 1778 PetscFunctionReturn(0); 1779 } 1780 /* -------------------------------------------------------------------------- */ 1781 1782 #undef __FUNCT__ 1783 #define __FUNCT__ "PCBDDCMatFETIDPGetRHS_BDDC" 1784 static PetscErrorCode PCBDDCMatFETIDPGetRHS_BDDC(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs) 1785 { 1786 FETIDPMat_ctx mat_ctx; 1787 Vec copy_standard_rhs; 1788 PC_IS* pcis; 1789 PC_BDDC* pcbddc; 1790 PetscErrorCode ierr; 1791 1792 PetscFunctionBegin; 1793 ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr); 1794 pcis = (PC_IS*)mat_ctx->pc->data; 1795 pcbddc = (PC_BDDC*)mat_ctx->pc->data; 1796 1797 /* 1798 change of basis for physical rhs if needed 1799 It also changes the rhs in case of dirichlet boundaries 1800 TODO: better management when FETIDP will have its own class 1801 */ 1802 ierr = VecDuplicate(standard_rhs,©_standard_rhs);CHKERRQ(ierr); 1803 ierr = VecCopy(standard_rhs,copy_standard_rhs);CHKERRQ(ierr); 1804 ierr = PCPreSolve_BDDC(mat_ctx->pc,NULL,copy_standard_rhs,NULL);CHKERRQ(ierr); 1805 /* store vectors for computation of fetidp final solution */ 1806 ierr = VecScatterBegin(pcis->global_to_D,copy_standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1807 ierr = VecScatterEnd(pcis->global_to_D,copy_standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1808 /* scale rhs since it should be unassembled */ 1809 /* TODO use counter scaling? (also below) */ 1810 ierr = VecScatterBegin(pcis->global_to_B,copy_standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1811 ierr = VecScatterEnd(pcis->global_to_B,copy_standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1812 /* Apply partition of unity */ 1813 ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr); 1814 /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,copy_standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */ 1815 if (!pcbddc->switch_static) { 1816 /* compute partially subassembled Schur complement right-hand side */ 1817 ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr); 1818 ierr = MatMult(pcis->A_BI,pcis->vec1_D,pcis->vec1_B);CHKERRQ(ierr); 1819 ierr = VecAXPY(mat_ctx->temp_solution_B,-1.0,pcis->vec1_B);CHKERRQ(ierr); 1820 ierr = VecSet(copy_standard_rhs,0.0);CHKERRQ(ierr); 1821 ierr = VecScatterBegin(pcis->global_to_B,mat_ctx->temp_solution_B,copy_standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1822 ierr = VecScatterEnd(pcis->global_to_B,mat_ctx->temp_solution_B,copy_standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1823 /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,copy_standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */ 1824 ierr = VecScatterBegin(pcis->global_to_B,copy_standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1825 ierr = VecScatterEnd(pcis->global_to_B,copy_standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1826 ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr); 1827 } 1828 ierr = VecDestroy(©_standard_rhs);CHKERRQ(ierr); 1829 /* BDDC rhs */ 1830 ierr = VecCopy(mat_ctx->temp_solution_B,pcis->vec1_B);CHKERRQ(ierr); 1831 if (pcbddc->switch_static) { 1832 ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr); 1833 } 1834 /* apply BDDC */ 1835 ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_FALSE);CHKERRQ(ierr); 1836 /* Application of B_delta and assembling of rhs for fetidp fluxes */ 1837 ierr = VecSet(fetidp_flux_rhs,0.0);CHKERRQ(ierr); 1838 ierr = MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local);CHKERRQ(ierr); 1839 ierr = VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1840 ierr = VecScatterEnd(mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1841 PetscFunctionReturn(0); 1842 } 1843 1844 #undef __FUNCT__ 1845 #define __FUNCT__ "PCBDDCMatFETIDPGetRHS" 1846 /*@ 1847 PCBDDCMatFETIDPGetRHS - Compute the right-hand side for FETI-DP linear system using the physical right-hand side 1848 1849 Collective 1850 1851 Input Parameters: 1852 + fetidp_mat - the FETI-DP matrix object obtained by a call to PCBDDCCreateFETIDPOperators 1853 - standard_rhs - the right-hand side of the original linear system 1854 1855 Output Parameters: 1856 . fetidp_flux_rhs - the right-hand side for the FETI-DP linear system 1857 1858 Level: developer 1859 1860 Notes: 1861 1862 .seealso: PCBDDC, PCBDDCCreateFETIDPOperators, PCBDDCMatFETIDPGetSolution 1863 @*/ 1864 PetscErrorCode PCBDDCMatFETIDPGetRHS(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs) 1865 { 1866 FETIDPMat_ctx mat_ctx; 1867 PetscErrorCode ierr; 1868 1869 PetscFunctionBegin; 1870 ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr); 1871 ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetRHS_C",(Mat,Vec,Vec),(fetidp_mat,standard_rhs,fetidp_flux_rhs));CHKERRQ(ierr); 1872 PetscFunctionReturn(0); 1873 } 1874 /* -------------------------------------------------------------------------- */ 1875 1876 #undef __FUNCT__ 1877 #define __FUNCT__ "PCBDDCMatFETIDPGetSolution_BDDC" 1878 static PetscErrorCode PCBDDCMatFETIDPGetSolution_BDDC(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol) 1879 { 1880 FETIDPMat_ctx mat_ctx; 1881 PC_IS* pcis; 1882 PC_BDDC* pcbddc; 1883 PetscErrorCode ierr; 1884 1885 PetscFunctionBegin; 1886 ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr); 1887 pcis = (PC_IS*)mat_ctx->pc->data; 1888 pcbddc = (PC_BDDC*)mat_ctx->pc->data; 1889 1890 /* apply B_delta^T */ 1891 ierr = VecScatterBegin(mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1892 ierr = VecScatterEnd (mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1893 ierr = MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr); 1894 /* compute rhs for BDDC application */ 1895 ierr = VecAYPX(pcis->vec1_B,-1.0,mat_ctx->temp_solution_B);CHKERRQ(ierr); 1896 if (pcbddc->switch_static) { 1897 ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr); 1898 } 1899 /* apply BDDC */ 1900 ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_FALSE);CHKERRQ(ierr); 1901 /* put values into standard global vector */ 1902 ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1903 ierr = VecScatterEnd (pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1904 if (!pcbddc->switch_static) { 1905 /* compute values into the interior if solved for the partially subassembled Schur complement */ 1906 ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec1_D);CHKERRQ(ierr); 1907 ierr = VecAXPY(mat_ctx->temp_solution_D,-1.0,pcis->vec1_D);CHKERRQ(ierr); 1908 ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr); 1909 } 1910 ierr = VecScatterBegin(pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1911 ierr = VecScatterEnd (pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1912 /* final change of basis if needed 1913 Is also sums the dirichlet part removed during RHS assembling */ 1914 ierr = PCPostSolve_BDDC(mat_ctx->pc,NULL,NULL,standard_sol);CHKERRQ(ierr); 1915 PetscFunctionReturn(0); 1916 } 1917 1918 #undef __FUNCT__ 1919 #define __FUNCT__ "PCBDDCMatFETIDPGetSolution" 1920 /*@ 1921 PCBDDCMatFETIDPGetSolution - Compute the physical solution using the solution of the FETI-DP linear system 1922 1923 Collective 1924 1925 Input Parameters: 1926 + fetidp_mat - the FETI-DP matrix obtained by a call to PCBDDCCreateFETIDPOperators 1927 - fetidp_flux_sol - the solution of the FETI-DP linear system 1928 1929 Output Parameters: 1930 . standard_sol - the solution defined on the physical domain 1931 1932 Level: developer 1933 1934 Notes: 1935 1936 .seealso: PCBDDC, PCBDDCCreateFETIDPOperators, PCBDDCMatFETIDPGetRHS 1937 @*/ 1938 PetscErrorCode PCBDDCMatFETIDPGetSolution(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol) 1939 { 1940 FETIDPMat_ctx mat_ctx; 1941 PetscErrorCode ierr; 1942 1943 PetscFunctionBegin; 1944 ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr); 1945 ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetSolution_C",(Mat,Vec,Vec),(fetidp_mat,fetidp_flux_sol,standard_sol));CHKERRQ(ierr); 1946 PetscFunctionReturn(0); 1947 } 1948 /* -------------------------------------------------------------------------- */ 1949 1950 extern PetscErrorCode FETIDPMatMult(Mat,Vec,Vec); 1951 extern PetscErrorCode FETIDPMatMultTranspose(Mat,Vec,Vec); 1952 extern PetscErrorCode PCBDDCDestroyFETIDPMat(Mat); 1953 extern PetscErrorCode FETIDPPCApply(PC,Vec,Vec); 1954 extern PetscErrorCode FETIDPPCApplyTranspose(PC,Vec,Vec); 1955 extern PetscErrorCode PCBDDCDestroyFETIDPPC(PC); 1956 1957 #undef __FUNCT__ 1958 #define __FUNCT__ "PCBDDCCreateFETIDPOperators_BDDC" 1959 static PetscErrorCode PCBDDCCreateFETIDPOperators_BDDC(PC pc, Mat *fetidp_mat, PC *fetidp_pc) 1960 { 1961 1962 FETIDPMat_ctx fetidpmat_ctx; 1963 Mat newmat; 1964 FETIDPPC_ctx fetidppc_ctx; 1965 PC newpc; 1966 MPI_Comm comm; 1967 PetscErrorCode ierr; 1968 1969 PetscFunctionBegin; 1970 ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 1971 /* FETIDP linear matrix */ 1972 ierr = PCBDDCCreateFETIDPMatContext(pc,&fetidpmat_ctx);CHKERRQ(ierr); 1973 ierr = PCBDDCSetupFETIDPMatContext(fetidpmat_ctx);CHKERRQ(ierr); 1974 ierr = MatCreateShell(comm,PETSC_DECIDE,PETSC_DECIDE,fetidpmat_ctx->n_lambda,fetidpmat_ctx->n_lambda,fetidpmat_ctx,&newmat);CHKERRQ(ierr); 1975 ierr = MatShellSetOperation(newmat,MATOP_MULT,(void (*)(void))FETIDPMatMult);CHKERRQ(ierr); 1976 ierr = MatShellSetOperation(newmat,MATOP_MULT_TRANSPOSE,(void (*)(void))FETIDPMatMultTranspose);CHKERRQ(ierr); 1977 ierr = MatShellSetOperation(newmat,MATOP_DESTROY,(void (*)(void))PCBDDCDestroyFETIDPMat);CHKERRQ(ierr); 1978 ierr = MatSetUp(newmat);CHKERRQ(ierr); 1979 /* FETIDP preconditioner */ 1980 ierr = PCBDDCCreateFETIDPPCContext(pc,&fetidppc_ctx);CHKERRQ(ierr); 1981 ierr = PCBDDCSetupFETIDPPCContext(newmat,fetidppc_ctx);CHKERRQ(ierr); 1982 ierr = PCCreate(comm,&newpc);CHKERRQ(ierr); 1983 ierr = PCSetType(newpc,PCSHELL);CHKERRQ(ierr); 1984 ierr = PCShellSetContext(newpc,fetidppc_ctx);CHKERRQ(ierr); 1985 ierr = PCShellSetApply(newpc,FETIDPPCApply);CHKERRQ(ierr); 1986 ierr = PCShellSetApplyTranspose(newpc,FETIDPPCApplyTranspose);CHKERRQ(ierr); 1987 ierr = PCShellSetDestroy(newpc,PCBDDCDestroyFETIDPPC);CHKERRQ(ierr); 1988 ierr = PCSetOperators(newpc,newmat,newmat);CHKERRQ(ierr); 1989 ierr = PCSetUp(newpc);CHKERRQ(ierr); 1990 /* return pointers for objects created */ 1991 *fetidp_mat=newmat; 1992 *fetidp_pc=newpc; 1993 PetscFunctionReturn(0); 1994 } 1995 1996 #undef __FUNCT__ 1997 #define __FUNCT__ "PCBDDCCreateFETIDPOperators" 1998 /*@ 1999 PCBDDCCreateFETIDPOperators - Create FETI-DP operators 2000 2001 Collective 2002 2003 Input Parameters: 2004 . pc - the BDDC preconditioning context (setup should have been called before) 2005 2006 Output Parameters: 2007 + fetidp_mat - shell FETI-DP matrix object 2008 - fetidp_pc - shell Dirichlet preconditioner for FETI-DP matrix 2009 2010 Options Database Keys: 2011 . -fetidp_fullyredundant <false> - use or not a fully redundant set of Lagrange multipliers 2012 2013 Level: developer 2014 2015 Notes: 2016 Currently the only operations provided for FETI-DP matrix are MatMult and MatMultTranspose 2017 2018 .seealso: PCBDDC, PCBDDCMatFETIDPGetRHS, PCBDDCMatFETIDPGetSolution 2019 @*/ 2020 PetscErrorCode PCBDDCCreateFETIDPOperators(PC pc, Mat *fetidp_mat, PC *fetidp_pc) 2021 { 2022 PetscErrorCode ierr; 2023 2024 PetscFunctionBegin; 2025 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 2026 if (pc->setupcalled) { 2027 ierr = PetscUseMethod(pc,"PCBDDCCreateFETIDPOperators_C",(PC,Mat*,PC*),(pc,fetidp_mat,fetidp_pc));CHKERRQ(ierr); 2028 } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"You must call PCSetup_BDDC() first \n"); 2029 PetscFunctionReturn(0); 2030 } 2031 /* -------------------------------------------------------------------------- */ 2032 /*MC 2033 PCBDDC - Balancing Domain Decomposition by Constraints. 2034 2035 An implementation of the BDDC preconditioner based on 2036 2037 .vb 2038 [1] C. R. Dohrmann. "An approximate BDDC preconditioner", Numerical Linear Algebra with Applications Volume 14, Issue 2, pages 149-168, March 2007 2039 [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 2040 [3] J. Mandel, B. Sousedik, C. R. Dohrmann. "Multispace and Multilevel BDDC", http://arxiv.org/abs/0712.3977 2041 [4] C. Pechstein and C. R. Dohrmann. "Modern domain decomposition methods BDDC, deluxe scaling, and an algebraic approach", Seminar talk, Linz, December 2013, http://people.ricam.oeaw.ac.at/c.pechstein/pechstein-bddc2013.pdf 2042 .ve 2043 2044 The matrix to be preconditioned (Pmat) must be of type MATIS. 2045 2046 Currently works with MATIS matrices with local matrices of type MATSEQAIJ, MATSEQBAIJ or MATSEQSBAIJ, either with real or complex numbers. 2047 2048 It also works with unsymmetric and indefinite problems. 2049 2050 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. 2051 2052 Approximate local solvers are automatically adapted for singular linear problems (see [1]) if the user has provided the nullspace using PCBDDCSetNullSpace() 2053 2054 Boundary nodes are split in vertices, edges and faces classes 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() 2055 Additional information on dofs can be provided by using PCBDDCSetDofsSplitting(), PCBDDCSetDirichletBoundaries(), PCBDDCSetNeumannBoundaries(), and PCBDDCSetPrimalVerticesLocalIS() 2056 2057 Constraints can be customized by attaching a MatNullSpace object to the MATIS matrix via MatSetNearNullSpace(). Non-singular modes are retained via SVD. 2058 2059 Change of basis is performed similarly to [2] when requested. When more than 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. 2060 User defined change of basis can be passed to PCBDDC by using PCBDDCSetChangeOfBasisMat() 2061 2062 The PETSc implementation also supports multilevel BDDC [3]. Coarse grids are partitioned using a MatPartitioning object. 2063 2064 Adaptive selection of primal constraints [4] is supported for SPD systems with high-contrast in the coefficients if MUMPS is present. Future versions of the code will also consider using MKL_PARDISO or PASTIX. 2065 2066 An experimental interface to the FETI-DP method is available. FETI-DP operators could be created using PCBDDCCreateFETIDPOperators(). A stand-alone class for the FETI-DP method will be provided in the next releases. 2067 Deluxe scaling is not supported yet for FETI-DP. 2068 2069 Options Database Keys (some of them, run with -h for a complete list): 2070 2071 . -pc_bddc_use_vertices <true> - use or not vertices in primal space 2072 . -pc_bddc_use_edges <true> - use or not edges in primal space 2073 . -pc_bddc_use_faces <false> - use or not faces in primal space 2074 . -pc_bddc_symmetric <true> - symmetric computation of primal basis functions. Specify false for unsymmetric problems 2075 . -pc_bddc_use_change_of_basis <false> - use change of basis approach (on edges only) 2076 . -pc_bddc_use_change_on_faces <false> - use change of basis approach on faces if change of basis has been requested 2077 . -pc_bddc_switch_static <false> - switches from M_2 (default) to M_3 operator (see reference article [1]) 2078 . -pc_bddc_levels <0> - maximum number of levels for multilevel 2079 . -pc_bddc_coarsening_ratio <8> - number of subdomains which will be aggregated together at the coarser level (e.g. H/h ratio at the coarser level, significative only in the multilevel case) 2080 . -pc_bddc_redistribute <0> - size of a subset of processors where the coarse problem will be remapped (the value is ignored if not at the coarsest level) 2081 . -pc_bddc_use_deluxe_scaling <false> - use deluxe scaling 2082 . -pc_bddc_schur_layers <-1> - select the economic version of deluxe scaling by specifying the number of layers (-1 corresponds to the original deluxe scaling) 2083 . -pc_bddc_adaptive_threshold <0.0> - when a value greater than one is specified, adaptive selection of constraints is performed on edges and faces (requires deluxe scaling and MUMPS installed) 2084 - -pc_bddc_check_level <0> - set verbosity level of debugging output 2085 2086 Options for Dirichlet, Neumann or coarse solver can be set with 2087 .vb 2088 -pc_bddc_dirichlet_ 2089 -pc_bddc_neumann_ 2090 -pc_bddc_coarse_ 2091 .ve 2092 e.g -pc_bddc_dirichlet_ksp_type richardson -pc_bddc_dirichlet_pc_type gamg. PCBDDC uses by default KPSPREONLY and PCLU. 2093 2094 When using a multilevel approach, solvers' options at the N-th level (N > 1) can be specified as 2095 .vb 2096 -pc_bddc_dirichlet_lN_ 2097 -pc_bddc_neumann_lN_ 2098 -pc_bddc_coarse_lN_ 2099 .ve 2100 Note that level number ranges from the finest (0) to the coarsest (N). 2101 In order to specify options for the BDDC operators at the coarser levels (and not for the solvers), prepend -pc_bddc_coarse_ or -pc_bddc_coarse_l to the option, e.g. 2102 .vb 2103 -pc_bddc_coarse_pc_bddc_adaptive_threshold 5 -pc_bddc_coarse_l1_pc_bddc_redistribute 3 2104 .ve 2105 will use a threshold of 5 for constraints' selection at the first coarse level and will redistribute the coarse problem of the first coarse level on 3 processors 2106 2107 Level: intermediate 2108 2109 Developer notes: 2110 2111 Contributed by Stefano Zampini 2112 2113 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, MATIS 2114 M*/ 2115 2116 #undef __FUNCT__ 2117 #define __FUNCT__ "PCCreate_BDDC" 2118 PETSC_EXTERN PetscErrorCode PCCreate_BDDC(PC pc) 2119 { 2120 PetscErrorCode ierr; 2121 PC_BDDC *pcbddc; 2122 2123 PetscFunctionBegin; 2124 /* Creates the private data structure for this preconditioner and attach it to the PC object. */ 2125 ierr = PetscNewLog(pc,&pcbddc);CHKERRQ(ierr); 2126 pc->data = (void*)pcbddc; 2127 2128 /* create PCIS data structure */ 2129 ierr = PCISCreate(pc);CHKERRQ(ierr); 2130 2131 /* BDDC customization */ 2132 pcbddc->use_local_adj = PETSC_TRUE; 2133 pcbddc->use_vertices = PETSC_TRUE; 2134 pcbddc->use_edges = PETSC_TRUE; 2135 pcbddc->use_faces = PETSC_FALSE; 2136 pcbddc->use_change_of_basis = PETSC_FALSE; 2137 pcbddc->use_change_on_faces = PETSC_FALSE; 2138 pcbddc->switch_static = PETSC_FALSE; 2139 pcbddc->use_nnsp_true = PETSC_FALSE; 2140 pcbddc->use_qr_single = PETSC_FALSE; 2141 pcbddc->symmetric_primal = PETSC_TRUE; 2142 pcbddc->benign_saddle_point = PETSC_FALSE; 2143 pcbddc->dbg_flag = 0; 2144 /* private */ 2145 pcbddc->local_primal_size = 0; 2146 pcbddc->local_primal_size_cc = 0; 2147 pcbddc->local_primal_ref_node = 0; 2148 pcbddc->local_primal_ref_mult = 0; 2149 pcbddc->n_vertices = 0; 2150 pcbddc->primal_indices_local_idxs = 0; 2151 pcbddc->recompute_topography = PETSC_FALSE; 2152 pcbddc->coarse_size = -1; 2153 pcbddc->new_primal_space = PETSC_FALSE; 2154 pcbddc->new_primal_space_local = PETSC_FALSE; 2155 pcbddc->global_primal_indices = 0; 2156 pcbddc->onearnullspace = 0; 2157 pcbddc->onearnullvecs_state = 0; 2158 pcbddc->user_primal_vertices = 0; 2159 pcbddc->NullSpace = 0; 2160 pcbddc->temp_solution = 0; 2161 pcbddc->original_rhs = 0; 2162 pcbddc->local_mat = 0; 2163 pcbddc->ChangeOfBasisMatrix = 0; 2164 pcbddc->user_ChangeOfBasisMatrix = 0; 2165 pcbddc->new_global_mat = 0; 2166 pcbddc->coarse_vec = 0; 2167 pcbddc->coarse_ksp = 0; 2168 pcbddc->coarse_phi_B = 0; 2169 pcbddc->coarse_phi_D = 0; 2170 pcbddc->coarse_psi_B = 0; 2171 pcbddc->coarse_psi_D = 0; 2172 pcbddc->vec1_P = 0; 2173 pcbddc->vec1_R = 0; 2174 pcbddc->vec2_R = 0; 2175 pcbddc->local_auxmat1 = 0; 2176 pcbddc->local_auxmat2 = 0; 2177 pcbddc->R_to_B = 0; 2178 pcbddc->R_to_D = 0; 2179 pcbddc->ksp_D = 0; 2180 pcbddc->ksp_R = 0; 2181 pcbddc->NeumannBoundaries = 0; 2182 pcbddc->NeumannBoundariesLocal = 0; 2183 pcbddc->DirichletBoundaries = 0; 2184 pcbddc->DirichletBoundariesLocal = 0; 2185 pcbddc->user_provided_isfordofs = PETSC_FALSE; 2186 pcbddc->n_ISForDofs = 0; 2187 pcbddc->n_ISForDofsLocal = 0; 2188 pcbddc->ISForDofs = 0; 2189 pcbddc->ISForDofsLocal = 0; 2190 pcbddc->ConstraintMatrix = 0; 2191 pcbddc->use_exact_dirichlet_trick = PETSC_TRUE; 2192 pcbddc->coarse_loc_to_glob = 0; 2193 pcbddc->coarsening_ratio = 8; 2194 pcbddc->coarse_adj_red = 0; 2195 pcbddc->current_level = 0; 2196 pcbddc->max_levels = 0; 2197 pcbddc->use_coarse_estimates = PETSC_FALSE; 2198 pcbddc->redistribute_coarse = 0; 2199 pcbddc->coarse_subassembling = 0; 2200 pcbddc->coarse_subassembling_init = 0; 2201 2202 /* benign subspace trick */ 2203 pcbddc->B0_ncol = 0; 2204 pcbddc->B0_cols = NULL; 2205 pcbddc->B0_vals = NULL; 2206 pcbddc->benign_change = 0; 2207 pcbddc->benign_vec = 0; 2208 pcbddc->benign_original_mat = 0; 2209 pcbddc->benign_sf = 0; 2210 pcbddc->benign_p0_lidx = -1; 2211 pcbddc->benign_p0_gidx = -1; 2212 2213 /* create local graph structure */ 2214 ierr = PCBDDCGraphCreate(&pcbddc->mat_graph);CHKERRQ(ierr); 2215 2216 /* scaling */ 2217 pcbddc->work_scaling = 0; 2218 pcbddc->use_deluxe_scaling = PETSC_FALSE; 2219 pcbddc->faster_deluxe = PETSC_FALSE; 2220 2221 /* create sub schurs structure */ 2222 ierr = PCBDDCSubSchursCreate(&pcbddc->sub_schurs);CHKERRQ(ierr); 2223 pcbddc->sub_schurs_rebuild = PETSC_FALSE; 2224 pcbddc->sub_schurs_layers = -1; 2225 pcbddc->sub_schurs_use_useradj = PETSC_FALSE; 2226 2227 pcbddc->computed_rowadj = PETSC_FALSE; 2228 2229 /* adaptivity */ 2230 pcbddc->adaptive_threshold = 0.0; 2231 pcbddc->adaptive_nmax = 0; 2232 pcbddc->adaptive_nmin = 0; 2233 2234 /* function pointers */ 2235 pc->ops->apply = PCApply_BDDC; 2236 pc->ops->applytranspose = PCApplyTranspose_BDDC; 2237 pc->ops->setup = PCSetUp_BDDC; 2238 pc->ops->destroy = PCDestroy_BDDC; 2239 pc->ops->setfromoptions = PCSetFromOptions_BDDC; 2240 pc->ops->view = 0; 2241 pc->ops->applyrichardson = 0; 2242 pc->ops->applysymmetricleft = 0; 2243 pc->ops->applysymmetricright = 0; 2244 pc->ops->presolve = PCPreSolve_BDDC; 2245 pc->ops->postsolve = PCPostSolve_BDDC; 2246 2247 /* composing function */ 2248 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetChangeOfBasisMat_C",PCBDDCSetChangeOfBasisMat_BDDC);CHKERRQ(ierr); 2249 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",PCBDDCSetPrimalVerticesLocalIS_BDDC);CHKERRQ(ierr); 2250 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",PCBDDCSetCoarseningRatio_BDDC);CHKERRQ(ierr); 2251 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",PCBDDCSetLevel_BDDC);CHKERRQ(ierr); 2252 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",PCBDDCSetUseExactDirichlet_BDDC);CHKERRQ(ierr); 2253 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",PCBDDCSetLevels_BDDC);CHKERRQ(ierr); 2254 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",PCBDDCSetNullSpace_BDDC);CHKERRQ(ierr); 2255 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",PCBDDCSetDirichletBoundaries_BDDC);CHKERRQ(ierr); 2256 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundariesLocal_C",PCBDDCSetDirichletBoundariesLocal_BDDC);CHKERRQ(ierr); 2257 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",PCBDDCSetNeumannBoundaries_BDDC);CHKERRQ(ierr); 2258 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundariesLocal_C",PCBDDCSetNeumannBoundariesLocal_BDDC);CHKERRQ(ierr); 2259 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",PCBDDCGetDirichletBoundaries_BDDC);CHKERRQ(ierr); 2260 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundariesLocal_C",PCBDDCGetDirichletBoundariesLocal_BDDC);CHKERRQ(ierr); 2261 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",PCBDDCGetNeumannBoundaries_BDDC);CHKERRQ(ierr); 2262 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundariesLocal_C",PCBDDCGetNeumannBoundariesLocal_BDDC);CHKERRQ(ierr); 2263 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",PCBDDCSetDofsSplitting_BDDC);CHKERRQ(ierr); 2264 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplittingLocal_C",PCBDDCSetDofsSplittingLocal_BDDC);CHKERRQ(ierr); 2265 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",PCBDDCSetLocalAdjacencyGraph_BDDC);CHKERRQ(ierr); 2266 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",PCBDDCCreateFETIDPOperators_BDDC);CHKERRQ(ierr); 2267 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",PCBDDCMatFETIDPGetRHS_BDDC);CHKERRQ(ierr); 2268 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",PCBDDCMatFETIDPGetSolution_BDDC);CHKERRQ(ierr); 2269 PetscFunctionReturn(0); 2270 } 2271 2272