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 PetscBool benign_correction_is_zero = PETSC_FALSE; 954 PetscBool iscg; 955 956 PetscFunctionBegin; 957 /* if we are working with cg, one dirichlet solve can be avoided during Krylov iterations */ 958 if (ksp) { 959 ierr = PetscObjectTypeCompare((PetscObject)ksp,KSPCG,&iscg);CHKERRQ(ierr); 960 if (!iscg) { 961 ierr = PCBDDCSetUseExactDirichlet(pc,PETSC_FALSE);CHKERRQ(ierr); 962 } 963 } 964 /* Creates parallel work vectors used in presolve */ 965 if (!pcbddc->original_rhs) { 966 ierr = VecDuplicate(pcis->vec1_global,&pcbddc->original_rhs);CHKERRQ(ierr); 967 } 968 if (!pcbddc->temp_solution) { 969 ierr = VecDuplicate(pcis->vec1_global,&pcbddc->temp_solution);CHKERRQ(ierr); 970 } 971 972 if (x) { 973 ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr); 974 used_vec = x; 975 } else { /* it can only happen when calling PCBDDCMatFETIDPGetRHS */ 976 ierr = PetscObjectReference((PetscObject)pcbddc->temp_solution);CHKERRQ(ierr); 977 used_vec = pcbddc->temp_solution; 978 ierr = VecSet(used_vec,0.0);CHKERRQ(ierr); 979 } 980 981 /* hack into ksp data structure since PCPreSolve comes earlier than setting to zero the guess in src/ksp/ksp/interface/itfunc.c */ 982 if (ksp) { 983 /* store the flag for the initial guess since it will be restored back during PCPostSolve_BDDC */ 984 ierr = KSPGetInitialGuessNonzero(ksp,&pcbddc->ksp_guess_nonzero);CHKERRQ(ierr); 985 if (!pcbddc->ksp_guess_nonzero) { 986 ierr = VecSet(used_vec,0.0);CHKERRQ(ierr); 987 } 988 } 989 990 pcbddc->rhs_change = PETSC_FALSE; 991 /* Take into account zeroed rows -> change rhs and store solution removed */ 992 if (rhs) { 993 IS dirIS = NULL; 994 995 /* DirichletBoundariesLocal may not be consistent among neighbours; gets a dirichlet dofs IS from graph (may be cached) */ 996 ierr = PCBDDCGraphGetDirichletDofs(pcbddc->mat_graph,&dirIS);CHKERRQ(ierr); 997 if (dirIS) { 998 Mat_IS *matis = (Mat_IS*)pc->pmat->data; 999 PetscInt dirsize,i,*is_indices; 1000 PetscScalar *array_x; 1001 const PetscScalar *array_diagonal; 1002 1003 ierr = MatGetDiagonal(pc->pmat,pcis->vec1_global);CHKERRQ(ierr); 1004 ierr = VecPointwiseDivide(pcis->vec1_global,rhs,pcis->vec1_global);CHKERRQ(ierr); 1005 ierr = VecScatterBegin(matis->rctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1006 ierr = VecScatterEnd(matis->rctx,pcis->vec1_global,pcis->vec2_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1007 ierr = VecScatterBegin(matis->rctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1008 ierr = VecScatterEnd(matis->rctx,used_vec,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1009 ierr = ISGetLocalSize(dirIS,&dirsize);CHKERRQ(ierr); 1010 ierr = VecGetArray(pcis->vec1_N,&array_x);CHKERRQ(ierr); 1011 ierr = VecGetArrayRead(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr); 1012 ierr = ISGetIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr); 1013 for (i=0; i<dirsize; i++) array_x[is_indices[i]] = array_diagonal[is_indices[i]]; 1014 ierr = ISRestoreIndices(dirIS,(const PetscInt**)&is_indices);CHKERRQ(ierr); 1015 ierr = VecRestoreArrayRead(pcis->vec2_N,&array_diagonal);CHKERRQ(ierr); 1016 ierr = VecRestoreArray(pcis->vec1_N,&array_x);CHKERRQ(ierr); 1017 ierr = VecScatterBegin(matis->rctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1018 ierr = VecScatterEnd(matis->rctx,pcis->vec1_N,used_vec,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1019 pcbddc->rhs_change = PETSC_TRUE; 1020 ierr = ISDestroy(&dirIS);CHKERRQ(ierr); 1021 } 1022 } 1023 1024 /* remove the computed solution or the initial guess from the rhs */ 1025 if (pcbddc->rhs_change || (ksp && pcbddc->ksp_guess_nonzero) ) { 1026 /* store the original rhs */ 1027 if (copy_rhs) { 1028 ierr = VecCopy(rhs,pcbddc->original_rhs);CHKERRQ(ierr); 1029 copy_rhs = PETSC_FALSE; 1030 } 1031 pcbddc->rhs_change = PETSC_TRUE; 1032 ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr); 1033 ierr = MatMultAdd(pc->mat,used_vec,rhs,rhs);CHKERRQ(ierr); 1034 ierr = VecScale(used_vec,-1.0);CHKERRQ(ierr); 1035 ierr = VecCopy(used_vec,pcbddc->temp_solution);CHKERRQ(ierr); 1036 if (ksp) { 1037 ierr = KSPSetInitialGuessNonzero(ksp,PETSC_FALSE);CHKERRQ(ierr); 1038 } 1039 } 1040 ierr = VecDestroy(&used_vec);CHKERRQ(ierr); 1041 1042 /* When using the benign trick: (TODO: what about FETI-DP?) 1043 - change rhs on pressures (iteration matrix is surely of type MATIS) 1044 - compute initial vector in benign space 1045 */ 1046 if (rhs && pcbddc->benign_saddle_point) { 1047 Mat_IS *matis = (Mat_IS*)(pc->mat->data); 1048 1049 /* store the original rhs */ 1050 if (copy_rhs) { 1051 ierr = VecCopy(rhs,pcbddc->original_rhs);CHKERRQ(ierr); 1052 copy_rhs = PETSC_FALSE; 1053 } 1054 1055 /* now change (locally) the basis */ 1056 ierr = VecScatterBegin(matis->rctx,rhs,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1057 ierr = VecScatterEnd(matis->rctx,rhs,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1058 if (pcbddc->benign_change) { 1059 ierr = MatMultTranspose(pcbddc->benign_change,pcis->vec1_N,pcis->vec2_N);CHKERRQ(ierr); 1060 ierr = VecScatterBegin(matis->rctx,pcis->vec2_N,rhs,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1061 ierr = VecScatterEnd(matis->rctx,pcis->vec2_N,rhs,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1062 /* change local iteration matrix */ 1063 ierr = MatPtAP(matis->A,pcbddc->benign_change,MAT_REUSE_MATRIX,2.0,&pcbddc->local_mat);CHKERRQ(ierr); 1064 ierr = MatDestroy(&pcbddc->benign_original_mat);CHKERRQ(ierr); 1065 ierr = PetscObjectReference((PetscObject)matis->A);CHKERRQ(ierr); 1066 pcbddc->benign_original_mat = matis->A; 1067 ierr = MatDestroy(&matis->A);CHKERRQ(ierr); 1068 ierr = PetscObjectReference((PetscObject)pcbddc->local_mat);CHKERRQ(ierr); 1069 matis->A = pcbddc->local_mat; 1070 ierr = MatDestroy(&pcbddc->local_mat);CHKERRQ(ierr); 1071 } else { 1072 ierr = VecScatterBegin(matis->rctx,pcis->vec1_N,rhs,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1073 ierr = VecScatterEnd(matis->rctx,pcis->vec1_N,rhs,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1074 } 1075 pcbddc->rhs_change = PETSC_TRUE; 1076 1077 /* compute u^*_h as in Xuemin Tu's thesis (see Section 4.8.1) */ 1078 /* TODO: what about Stokes? */ 1079 if (!pcbddc->benign_vec) { 1080 ierr = VecDuplicate(rhs,&pcbddc->benign_vec);CHKERRQ(ierr); 1081 } 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 if (pcbddc->benign_null && iscg) { /* this is a workaround, need to understand more */ 1091 PetscBool iszero_l; 1092 1093 iszero_l = PetscAbsScalar(pcbddc->benign_p0) < PETSC_SMALL ? PETSC_TRUE : PETSC_FALSE; 1094 ierr = MPI_Allreduce(&iszero_l,&benign_correction_is_zero,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr); 1095 } 1096 if (!benign_correction_is_zero) { 1097 ierr = VecSet(pcis->vec1_global,0.);CHKERRQ(ierr); 1098 ierr = PCBDDCBenignGetOrSetP0(pc,pcis->vec1_global,PETSC_FALSE);CHKERRQ(ierr); 1099 ierr = PCApply_BDDC(pc,pcis->vec1_global,pcbddc->benign_vec);CHKERRQ(ierr); 1100 pcbddc->benign_p0 = 0.; 1101 ierr = PCBDDCBenignGetOrSetP0(pc,pcbddc->benign_vec,PETSC_FALSE);CHKERRQ(ierr); 1102 } 1103 } 1104 1105 /* change rhs and iteration matrix if using the change of basis */ 1106 if (pcbddc->ChangeOfBasisMatrix) { 1107 PCBDDCChange_ctx change_ctx; 1108 1109 /* get change ctx */ 1110 ierr = MatShellGetContext(pcbddc->new_global_mat,&change_ctx);CHKERRQ(ierr); 1111 1112 /* set current iteration matrix inside change context (change of basis has been already set into the ctx during PCSetUp) */ 1113 ierr = MatDestroy(&change_ctx->original_mat);CHKERRQ(ierr); 1114 ierr = PetscObjectReference((PetscObject)pc->mat);CHKERRQ(ierr); 1115 change_ctx->original_mat = pc->mat; 1116 1117 /* change iteration matrix */ 1118 ierr = MatDestroy(&pc->mat);CHKERRQ(ierr); 1119 ierr = PetscObjectReference((PetscObject)pcbddc->new_global_mat);CHKERRQ(ierr); 1120 pc->mat = pcbddc->new_global_mat; 1121 1122 /* store the original rhs */ 1123 if (copy_rhs) { 1124 ierr = VecCopy(rhs,pcbddc->original_rhs);CHKERRQ(ierr); 1125 copy_rhs = PETSC_FALSE; 1126 } 1127 1128 /* change rhs */ 1129 ierr = MatMultTranspose(change_ctx->global_change,rhs,pcis->vec1_global);CHKERRQ(ierr); 1130 ierr = VecCopy(pcis->vec1_global,rhs);CHKERRQ(ierr); 1131 pcbddc->rhs_change = PETSC_TRUE; 1132 } 1133 1134 /* remove non-benign solution from the rhs */ 1135 if (pcbddc->benign_saddle_point) { 1136 /* store the original rhs */ 1137 if (copy_rhs) { 1138 ierr = VecCopy(rhs,pcbddc->original_rhs);CHKERRQ(ierr); 1139 copy_rhs = PETSC_FALSE; 1140 } 1141 if (benign_correction_is_zero) { /* still need to understand why it works great */ 1142 ierr = PCBDDCSetUseExactDirichlet(pc,PETSC_FALSE);CHKERRQ(ierr); 1143 ierr = PCApply_BDDC(pc,rhs,pcbddc->benign_vec);CHKERRQ(ierr); 1144 } 1145 ierr = VecScale(pcbddc->benign_vec,-1.0);CHKERRQ(ierr); 1146 ierr = MatMultAdd(pc->mat,pcbddc->benign_vec,rhs,rhs);CHKERRQ(ierr); 1147 pcbddc->rhs_change = PETSC_TRUE; 1148 } 1149 1150 /* set initial guess if using PCG */ 1151 if (x && pcbddc->use_exact_dirichlet_trick) { 1152 ierr = VecSet(x,0.0);CHKERRQ(ierr); 1153 ierr = VecScatterBegin(pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1154 ierr = VecScatterEnd(pcis->global_to_D,rhs,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1155 ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr); 1156 ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1157 ierr = VecScatterEnd(pcis->global_to_D,pcis->vec2_D,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1158 if (ksp) { 1159 ierr = KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);CHKERRQ(ierr); 1160 } 1161 } 1162 if (pcbddc->benign_null) { 1163 MatNullSpace null_space; 1164 Vec nullv; 1165 PetscBool isnull; 1166 1167 pcbddc->benign_p0 = 1.; 1168 ierr = VecDuplicate(pcis->vec1_global,&nullv);CHKERRQ(ierr); 1169 ierr = VecSet(nullv,0.);CHKERRQ(ierr); 1170 ierr = PCBDDCBenignGetOrSetP0(pc,nullv,PETSC_FALSE);CHKERRQ(ierr); 1171 ierr = VecNormalize(nullv,NULL);CHKERRQ(ierr); 1172 ierr = MatNullSpaceCreate(PetscObjectComm((PetscObject)pc),PETSC_FALSE,1,&nullv,&null_space);CHKERRQ(ierr); 1173 ierr = MatNullSpaceTest(null_space,pc->mat,&isnull);CHKERRQ(ierr); 1174 ierr = MatSetNullSpace(pc->mat,null_space);CHKERRQ(ierr); 1175 ierr = MatNullSpaceDestroy(&null_space);CHKERRQ(ierr); 1176 ierr = VecDestroy(&nullv);CHKERRQ(ierr); 1177 } 1178 1179 1180 /* remove nullspace if present */ 1181 if (ksp && x && pcbddc->NullSpace) { 1182 ierr = MatNullSpaceRemove(pcbddc->NullSpace,x);CHKERRQ(ierr); 1183 /* store the original rhs */ 1184 if (copy_rhs) { 1185 ierr = VecCopy(rhs,pcbddc->original_rhs);CHKERRQ(ierr); 1186 copy_rhs = PETSC_FALSE; 1187 } 1188 pcbddc->rhs_change = PETSC_TRUE; 1189 ierr = MatNullSpaceRemove(pcbddc->NullSpace,rhs);CHKERRQ(ierr); 1190 } 1191 PetscFunctionReturn(0); 1192 } 1193 1194 /* -------------------------------------------------------------------------- */ 1195 #undef __FUNCT__ 1196 #define __FUNCT__ "PCPostSolve_BDDC" 1197 /* -------------------------------------------------------------------------- */ 1198 /* 1199 PCPostSolve_BDDC - Changes the computed solution if a transformation of basis 1200 approach has been selected. Also, restores rhs to its original state. 1201 1202 Input Parameter: 1203 + pc - the preconditioner contex 1204 1205 Application Interface Routine: PCPostSolve() 1206 1207 Notes: 1208 The interface routine PCPostSolve() is not usually called directly by 1209 the user, but instead is called by KSPSolve(). 1210 */ 1211 static PetscErrorCode PCPostSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x) 1212 { 1213 PetscErrorCode ierr; 1214 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 1215 1216 PetscFunctionBegin; 1217 if (pcbddc->ChangeOfBasisMatrix) { 1218 PCBDDCChange_ctx change_ctx; 1219 1220 /* get change ctx */ 1221 ierr = MatShellGetContext(pcbddc->new_global_mat,&change_ctx);CHKERRQ(ierr); 1222 1223 /* restore iteration matrix */ 1224 ierr = MatDestroy(&pc->mat);CHKERRQ(ierr); 1225 ierr = PetscObjectReference((PetscObject)change_ctx->original_mat);CHKERRQ(ierr); 1226 pc->mat = change_ctx->original_mat; 1227 } 1228 1229 /* need to restore the local matrices */ 1230 if (pcbddc->benign_change) { 1231 Mat_IS *matis = (Mat_IS*)pc->mat->data; 1232 1233 ierr = MatDestroy(&matis->A);CHKERRQ(ierr); 1234 ierr = PetscObjectReference((PetscObject)pcbddc->benign_original_mat);CHKERRQ(ierr); 1235 matis->A = pcbddc->benign_original_mat; 1236 ierr = MatDestroy(&pcbddc->benign_original_mat);CHKERRQ(ierr); 1237 } 1238 1239 /* get solution in original basis */ 1240 if (x) { 1241 PC_IS *pcis = (PC_IS*)(pc->data); 1242 1243 /* restore solution on pressures */ 1244 if (pcbddc->benign_saddle_point) { 1245 Mat_IS *matis = (Mat_IS*)pc->mat->data; 1246 1247 /* add non-benign solution */ 1248 ierr = VecAXPY(x,-1.0,pcbddc->benign_vec);CHKERRQ(ierr); 1249 1250 /* change basis on pressures for x */ 1251 ierr = VecScatterBegin(matis->rctx,x,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1252 ierr = VecScatterEnd(matis->rctx,x,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1253 if (pcbddc->benign_change) { 1254 ierr = MatMult(pcbddc->benign_change,pcis->vec1_N,pcis->vec2_N);CHKERRQ(ierr); 1255 ierr = VecScatterBegin(matis->rctx,pcis->vec2_N,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1256 ierr = VecScatterEnd(matis->rctx,pcis->vec2_N,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1257 } else { 1258 ierr = VecScatterBegin(matis->rctx,pcis->vec1_N,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1259 ierr = VecScatterEnd(matis->rctx,pcis->vec1_N,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1260 } 1261 } 1262 1263 /* change basis on x */ 1264 if (pcbddc->ChangeOfBasisMatrix) { 1265 PCBDDCChange_ctx change_ctx; 1266 1267 ierr = MatShellGetContext(pcbddc->new_global_mat,&change_ctx);CHKERRQ(ierr); 1268 ierr = MatMult(change_ctx->global_change,x,pcis->vec1_global);CHKERRQ(ierr); 1269 ierr = VecCopy(pcis->vec1_global,x);CHKERRQ(ierr); 1270 } 1271 } 1272 1273 /* add solution removed in presolve */ 1274 if (x && pcbddc->rhs_change) { 1275 ierr = VecAXPY(x,1.0,pcbddc->temp_solution);CHKERRQ(ierr); 1276 } 1277 1278 /* restore rhs to its original state */ 1279 if (rhs && pcbddc->rhs_change) { 1280 ierr = VecCopy(pcbddc->original_rhs,rhs);CHKERRQ(ierr); 1281 } 1282 pcbddc->rhs_change = PETSC_FALSE; 1283 1284 /* restore ksp guess state */ 1285 if (ksp) { 1286 ierr = KSPSetInitialGuessNonzero(ksp,pcbddc->ksp_guess_nonzero);CHKERRQ(ierr); 1287 } 1288 PetscFunctionReturn(0); 1289 } 1290 /* -------------------------------------------------------------------------- */ 1291 #undef __FUNCT__ 1292 #define __FUNCT__ "PCSetUp_BDDC" 1293 /* -------------------------------------------------------------------------- */ 1294 /* 1295 PCSetUp_BDDC - Prepares for the use of the BDDC preconditioner 1296 by setting data structures and options. 1297 1298 Input Parameter: 1299 + pc - the preconditioner context 1300 1301 Application Interface Routine: PCSetUp() 1302 1303 Notes: 1304 The interface routine PCSetUp() is not usually called directly by 1305 the user, but instead is called by PCApply() if necessary. 1306 */ 1307 PetscErrorCode PCSetUp_BDDC(PC pc) 1308 { 1309 PetscErrorCode ierr; 1310 PC_BDDC* pcbddc = (PC_BDDC*)pc->data; 1311 Mat_IS* matis; 1312 MatNullSpace nearnullspace; 1313 IS zerodiag = NULL; 1314 PetscInt nrows,ncols; 1315 PetscBool computetopography,computesolvers,computesubschurs; 1316 PetscBool computeconstraintsmatrix; 1317 PetscBool new_nearnullspace_provided,ismatis; 1318 1319 PetscFunctionBegin; 1320 ierr = PetscObjectTypeCompare((PetscObject)pc->pmat,MATIS,&ismatis);CHKERRQ(ierr); 1321 if (!ismatis) { 1322 SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"PCBDDC preconditioner requires matrix of type MATIS"); 1323 } 1324 ierr = MatGetSize(pc->pmat,&nrows,&ncols);CHKERRQ(ierr); 1325 if (nrows != ncols) { 1326 SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"PCBDDC preconditioner requires a square preconditioning matrix"); 1327 } 1328 matis = (Mat_IS*)pc->pmat->data; 1329 /* the following lines of code should be replaced by a better logic between PCIS, PCNN, PCBDDC and other future nonoverlapping preconditioners */ 1330 /* For BDDC we need to define a local "Neumann" problem different to that defined in PCISSetup 1331 Also, BDDC directly build the Dirichlet problem */ 1332 /* split work */ 1333 if (pc->setupcalled) { 1334 if (pc->flag == SAME_NONZERO_PATTERN) { 1335 computetopography = PETSC_FALSE; 1336 computesolvers = PETSC_TRUE; 1337 } else { /* DIFFERENT_NONZERO_PATTERN */ 1338 computetopography = PETSC_TRUE; 1339 computesolvers = PETSC_TRUE; 1340 } 1341 } else { 1342 computetopography = PETSC_TRUE; 1343 computesolvers = PETSC_TRUE; 1344 } 1345 if (pcbddc->recompute_topography) { 1346 computetopography = PETSC_TRUE; 1347 } 1348 computeconstraintsmatrix = PETSC_FALSE; 1349 1350 /* check parameters' compatibility */ 1351 if (pcbddc->adaptive_threshold > 0.0 && !pcbddc->use_deluxe_scaling) { 1352 SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Cannot compute adaptive constraints without deluxe scaling. Rerun with -pc_bddc_use_deluxe_scaling"); 1353 } 1354 pcbddc->adaptive_selection = (PetscBool)(pcbddc->adaptive_threshold > 0.0 && pcbddc->use_deluxe_scaling); 1355 if (pcbddc->adaptive_selection) pcbddc->use_faces = PETSC_TRUE; 1356 1357 computesubschurs = (PetscBool)(pcbddc->adaptive_selection || pcbddc->use_deluxe_scaling); 1358 if (pcbddc->faster_deluxe && pcbddc->adaptive_selection && pcbddc->use_change_of_basis) { 1359 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"); 1360 } 1361 1362 /* check if the iteration matrix is of type MATIS in case the benign trick has been requested */ 1363 ierr = PetscObjectTypeCompare((PetscObject)pc->mat,MATIS,&ismatis);CHKERRQ(ierr); 1364 if (pcbddc->benign_saddle_point && !ismatis) { 1365 SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"PCBDDC preconditioner with benign subspace trick requires the iteration matrix to be of type MATIS"); 1366 } 1367 1368 /* Get stdout for dbg */ 1369 if (pcbddc->dbg_flag) { 1370 if (!pcbddc->dbg_viewer) { 1371 pcbddc->dbg_viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)pc)); 1372 ierr = PetscViewerASCIISynchronizedAllow(pcbddc->dbg_viewer,PETSC_TRUE);CHKERRQ(ierr); 1373 } 1374 ierr = PetscViewerASCIIAddTab(pcbddc->dbg_viewer,2*pcbddc->current_level);CHKERRQ(ierr); 1375 } 1376 1377 if (pcbddc->user_ChangeOfBasisMatrix) { 1378 /* use_change_of_basis flag is used to automatically compute a change of basis from constraints */ 1379 pcbddc->use_change_of_basis = PETSC_FALSE; 1380 ierr = PCBDDCComputeLocalMatrix(pc,pcbddc->user_ChangeOfBasisMatrix);CHKERRQ(ierr); 1381 } else { 1382 ierr = MatDestroy(&pcbddc->local_mat);CHKERRQ(ierr); 1383 ierr = PetscObjectReference((PetscObject)matis->A);CHKERRQ(ierr); 1384 pcbddc->local_mat = matis->A; 1385 } 1386 1387 /* change basis on local pressures (aka zerodiag dofs) */ 1388 if (pcbddc->benign_saddle_point) { 1389 /* detect local saddle point and change the basis in pcbddc->local_mat */ 1390 ierr = PCBDDCBenignDetectSaddlePoint(pc,&zerodiag);CHKERRQ(ierr); 1391 /* pop B0 mat from pcbddc->local_mat */ 1392 ierr = PCBDDCBenignPopOrPushB0(pc,PETSC_TRUE);CHKERRQ(ierr); 1393 } 1394 1395 /* propagate relevant information */ 1396 /* workaround for reals */ 1397 #if !defined(PETSC_USE_COMPLEX) 1398 if (matis->A->symmetric_set) { 1399 ierr = MatSetOption(pcbddc->local_mat,MAT_HERMITIAN,matis->A->symmetric);CHKERRQ(ierr); 1400 } 1401 #endif 1402 if (matis->A->symmetric_set) { 1403 ierr = MatSetOption(pcbddc->local_mat,MAT_SYMMETRIC,matis->A->symmetric);CHKERRQ(ierr); 1404 } 1405 if (matis->A->spd_set) { 1406 ierr = MatSetOption(pcbddc->local_mat,MAT_SPD,matis->A->spd);CHKERRQ(ierr); 1407 } 1408 1409 /* Set up all the "iterative substructuring" common block without computing solvers */ 1410 { 1411 Mat temp_mat; 1412 1413 temp_mat = matis->A; 1414 matis->A = pcbddc->local_mat; 1415 ierr = PCISSetUp(pc,PETSC_FALSE);CHKERRQ(ierr); 1416 pcbddc->local_mat = matis->A; 1417 matis->A = temp_mat; 1418 } 1419 1420 /* Analyze interface */ 1421 if (computetopography) { 1422 ierr = PCBDDCAnalyzeInterface(pc);CHKERRQ(ierr); 1423 computeconstraintsmatrix = PETSC_TRUE; 1424 } 1425 1426 /* check existence of a divergence free extension, i.e. 1427 b(v_I,p_0) = 0 for all v_I (raise error if not). 1428 Also, check that PCBDDCBenignGetOrSetP0 works */ 1429 #if defined(PETSC_USE_DEBUG) 1430 if (pcbddc->benign_saddle_point) { 1431 ierr = PCBDDCBenignCheck(pc,zerodiag);CHKERRQ(ierr); 1432 } 1433 #endif 1434 1435 /* Setup local dirichlet solver ksp_D and sub_schurs solvers */ 1436 if (computesolvers) { 1437 PCBDDCSubSchurs sub_schurs=pcbddc->sub_schurs; 1438 1439 if (computesubschurs && computetopography) { 1440 ierr = PCBDDCInitSubSchurs(pc);CHKERRQ(ierr); 1441 } 1442 if (sub_schurs->use_mumps) { 1443 if (computesubschurs) { 1444 ierr = PCBDDCSetUpSubSchurs(pc);CHKERRQ(ierr); 1445 } 1446 ierr = PCBDDCSetUpLocalSolvers(pc,PETSC_TRUE,PETSC_FALSE);CHKERRQ(ierr); 1447 } else { 1448 ierr = PCBDDCSetUpLocalSolvers(pc,PETSC_TRUE,PETSC_FALSE);CHKERRQ(ierr); 1449 if (computesubschurs) { 1450 ierr = PCBDDCSetUpSubSchurs(pc);CHKERRQ(ierr); 1451 } 1452 } 1453 if (pcbddc->adaptive_selection) { 1454 ierr = PCBDDCAdaptiveSelection(pc);CHKERRQ(ierr); 1455 computeconstraintsmatrix = PETSC_TRUE; 1456 } 1457 } 1458 1459 /* infer if NullSpace object attached to Mat via MatSetNearNullSpace has changed */ 1460 new_nearnullspace_provided = PETSC_FALSE; 1461 ierr = MatGetNearNullSpace(pc->pmat,&nearnullspace);CHKERRQ(ierr); 1462 if (pcbddc->onearnullspace) { /* already used nearnullspace */ 1463 if (!nearnullspace) { /* near null space attached to mat has been destroyed */ 1464 new_nearnullspace_provided = PETSC_TRUE; 1465 } else { 1466 /* determine if the two nullspaces are different (should be lightweight) */ 1467 if (nearnullspace != pcbddc->onearnullspace) { 1468 new_nearnullspace_provided = PETSC_TRUE; 1469 } else { /* maybe the user has changed the content of the nearnullspace so check vectors ObjectStateId */ 1470 PetscInt i; 1471 const Vec *nearnullvecs; 1472 PetscObjectState state; 1473 PetscInt nnsp_size; 1474 ierr = MatNullSpaceGetVecs(nearnullspace,NULL,&nnsp_size,&nearnullvecs);CHKERRQ(ierr); 1475 for (i=0;i<nnsp_size;i++) { 1476 ierr = PetscObjectStateGet((PetscObject)nearnullvecs[i],&state);CHKERRQ(ierr); 1477 if (pcbddc->onearnullvecs_state[i] != state) { 1478 new_nearnullspace_provided = PETSC_TRUE; 1479 break; 1480 } 1481 } 1482 } 1483 } 1484 } else { 1485 if (!nearnullspace) { /* both nearnullspaces are null */ 1486 new_nearnullspace_provided = PETSC_FALSE; 1487 } else { /* nearnullspace attached later */ 1488 new_nearnullspace_provided = PETSC_TRUE; 1489 } 1490 } 1491 1492 /* Setup constraints and related work vectors */ 1493 /* reset primal space flags */ 1494 pcbddc->new_primal_space = PETSC_FALSE; 1495 pcbddc->new_primal_space_local = PETSC_FALSE; 1496 if (computeconstraintsmatrix || new_nearnullspace_provided) { 1497 /* It also sets the primal space flags */ 1498 ierr = PCBDDCConstraintsSetUp(pc);CHKERRQ(ierr); 1499 /* Allocate needed local vectors (which depends on quantities defined during ConstraintsSetUp) */ 1500 ierr = PCBDDCSetUpLocalWorkVectors(pc);CHKERRQ(ierr); 1501 } 1502 1503 if (computesolvers || pcbddc->new_primal_space) { 1504 if (pcbddc->use_change_of_basis) { 1505 PC_IS *pcis = (PC_IS*)(pc->data); 1506 Mat temp_mat = NULL; 1507 1508 if (zerodiag) { 1509 /* insert B0 in pcbddc->local_mat */ 1510 ierr = PCBDDCBenignPopOrPushB0(pc,PETSC_FALSE);CHKERRQ(ierr); 1511 /* swap local matrices */ 1512 ierr = MatISGetLocalMat(pc->pmat,&temp_mat);CHKERRQ(ierr); 1513 ierr = PetscObjectReference((PetscObject)temp_mat);CHKERRQ(ierr); 1514 ierr = MatISSetLocalMat(pc->pmat,pcbddc->local_mat);CHKERRQ(ierr); 1515 ierr = MatDestroy(&pcbddc->local_mat);CHKERRQ(ierr); 1516 } 1517 ierr = PCBDDCComputeLocalMatrix(pc,pcbddc->ChangeOfBasisMatrix);CHKERRQ(ierr); 1518 if (zerodiag) { 1519 /* restore original matrix */ 1520 ierr = MatISSetLocalMat(pc->pmat,temp_mat);CHKERRQ(ierr); 1521 ierr = PetscObjectDereference((PetscObject)temp_mat);CHKERRQ(ierr); 1522 /* pop B0 from pcbddc->local_mat */ 1523 ierr = PCBDDCBenignPopOrPushB0(pc,PETSC_TRUE);CHKERRQ(ierr); 1524 } 1525 /* get submatrices */ 1526 ierr = MatDestroy(&pcis->A_IB);CHKERRQ(ierr); 1527 ierr = MatDestroy(&pcis->A_BI);CHKERRQ(ierr); 1528 ierr = MatDestroy(&pcis->A_BB);CHKERRQ(ierr); 1529 ierr = MatGetSubMatrix(pcbddc->local_mat,pcis->is_B_local,pcis->is_B_local,MAT_INITIAL_MATRIX,&pcis->A_BB);CHKERRQ(ierr); 1530 ierr = MatGetSubMatrix(pcbddc->local_mat,pcis->is_I_local,pcis->is_B_local,MAT_INITIAL_MATRIX,&pcis->A_IB);CHKERRQ(ierr); 1531 ierr = MatGetSubMatrix(pcbddc->local_mat,pcis->is_B_local,pcis->is_I_local,MAT_INITIAL_MATRIX,&pcis->A_BI);CHKERRQ(ierr); 1532 /* set flag in pcis to not reuse submatrices during PCISCreate */ 1533 pcis->reusesubmatrices = PETSC_FALSE; 1534 } else if (!pcbddc->user_ChangeOfBasisMatrix && !pcbddc->benign_change) { 1535 ierr = MatDestroy(&pcbddc->local_mat);CHKERRQ(ierr); 1536 ierr = PetscObjectReference((PetscObject)matis->A);CHKERRQ(ierr); 1537 pcbddc->local_mat = matis->A; 1538 } 1539 /* SetUp coarse and local Neumann solvers */ 1540 ierr = PCBDDCSetUpSolvers(pc);CHKERRQ(ierr); 1541 /* SetUp Scaling operator */ 1542 ierr = PCBDDCScalingSetUp(pc);CHKERRQ(ierr); 1543 } 1544 ierr = ISDestroy(&zerodiag);CHKERRQ(ierr); 1545 1546 if (pcbddc->dbg_flag) { 1547 ierr = PetscViewerASCIISubtractTab(pcbddc->dbg_viewer,2*pcbddc->current_level);CHKERRQ(ierr); 1548 } 1549 PetscFunctionReturn(0); 1550 } 1551 1552 /* -------------------------------------------------------------------------- */ 1553 /* 1554 PCApply_BDDC - Applies the BDDC operator to a vector. 1555 1556 Input Parameters: 1557 + pc - the preconditioner context 1558 - r - input vector (global) 1559 1560 Output Parameter: 1561 . z - output vector (global) 1562 1563 Application Interface Routine: PCApply() 1564 */ 1565 #undef __FUNCT__ 1566 #define __FUNCT__ "PCApply_BDDC" 1567 PetscErrorCode PCApply_BDDC(PC pc,Vec r,Vec z) 1568 { 1569 PC_IS *pcis = (PC_IS*)(pc->data); 1570 PC_BDDC *pcbddc = (PC_BDDC*)(pc->data); 1571 PetscInt n_B = pcis->n_B, n_D = pcis->n - n_B; 1572 PetscErrorCode ierr; 1573 const PetscScalar one = 1.0; 1574 const PetscScalar m_one = -1.0; 1575 const PetscScalar zero = 0.0; 1576 1577 /* This code is similar to that provided in nn.c for PCNN 1578 NN interface preconditioner changed to BDDC 1579 Added support for M_3 preconditioner in the reference article (code is active if pcbddc->switch_static == PETSC_TRUE) */ 1580 1581 PetscFunctionBegin; 1582 if (pcbddc->benign_saddle_point) { /* get p0 from r */ 1583 ierr = PCBDDCBenignGetOrSetP0(pc,r,PETSC_TRUE);CHKERRQ(ierr); 1584 } 1585 if (!pcbddc->use_exact_dirichlet_trick) { 1586 ierr = VecCopy(r,z);CHKERRQ(ierr); 1587 /* First Dirichlet solve */ 1588 ierr = VecScatterBegin(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1589 ierr = VecScatterEnd(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1590 /* 1591 Assembling right hand side for BDDC operator 1592 - pcis->vec1_D for the Dirichlet part (if needed, i.e. pcbddc->switch_static == PETSC_TRUE) 1593 - pcis->vec1_B the interface part of the global vector z 1594 */ 1595 if (n_D) { 1596 ierr = KSPSolve(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr); 1597 ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr); 1598 if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec2_D,pcis->vec1_D,pcis->vec1_D);CHKERRQ(ierr); } 1599 ierr = MatMult(pcis->A_BI,pcis->vec2_D,pcis->vec1_B);CHKERRQ(ierr); 1600 } else { 1601 ierr = VecSet(pcis->vec1_B,zero);CHKERRQ(ierr); 1602 } 1603 ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1604 ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1605 ierr = PCBDDCScalingRestriction(pc,z,pcis->vec1_B);CHKERRQ(ierr); 1606 } else { 1607 if (pcbddc->switch_static) { 1608 ierr = VecSet(pcis->vec1_D,zero);CHKERRQ(ierr); 1609 } 1610 ierr = PCBDDCScalingRestriction(pc,r,pcis->vec1_B);CHKERRQ(ierr); 1611 } 1612 1613 /* Apply interface preconditioner 1614 input/output vecs: pcis->vec1_B and pcis->vec1_D */ 1615 ierr = PCBDDCApplyInterfacePreconditioner(pc,PETSC_FALSE);CHKERRQ(ierr); 1616 1617 /* Apply transpose of partition of unity operator */ 1618 ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,z);CHKERRQ(ierr); 1619 1620 /* Second Dirichlet solve and assembling of output */ 1621 ierr = VecScatterBegin(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1622 ierr = VecScatterEnd(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1623 if (n_B) { 1624 ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec3_D);CHKERRQ(ierr); 1625 if (pcbddc->switch_static) { ierr = MatMultAdd(pcis->A_II,pcis->vec1_D,pcis->vec3_D,pcis->vec3_D);CHKERRQ(ierr); } 1626 } else if (pcbddc->switch_static) { 1627 ierr = MatMult(pcis->A_II,pcis->vec1_D,pcis->vec3_D);CHKERRQ(ierr); 1628 } 1629 ierr = KSPSolve(pcbddc->ksp_D,pcis->vec3_D,pcis->vec4_D);CHKERRQ(ierr); 1630 1631 if (!pcbddc->use_exact_dirichlet_trick) { 1632 if (pcbddc->switch_static) { 1633 ierr = VecAXPBYPCZ(pcis->vec2_D,m_one,one,m_one,pcis->vec4_D,pcis->vec1_D);CHKERRQ(ierr); 1634 } else { 1635 ierr = VecAXPBY(pcis->vec2_D,m_one,m_one,pcis->vec4_D);CHKERRQ(ierr); 1636 } 1637 ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1638 ierr = VecScatterEnd(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1639 } else { 1640 if (pcbddc->switch_static) { 1641 ierr = VecAXPBY(pcis->vec4_D,one,m_one,pcis->vec1_D);CHKERRQ(ierr); 1642 } else { 1643 ierr = VecScale(pcis->vec4_D,m_one);CHKERRQ(ierr); 1644 } 1645 ierr = VecScatterBegin(pcis->global_to_D,pcis->vec4_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1646 ierr = VecScatterEnd(pcis->global_to_D,pcis->vec4_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1647 } 1648 1649 if (pcbddc->benign_saddle_point) { /* set p0 (computed in PCBDDCApplyInterface) */ 1650 ierr = PCBDDCBenignGetOrSetP0(pc,z,PETSC_FALSE);CHKERRQ(ierr); 1651 } 1652 PetscFunctionReturn(0); 1653 } 1654 1655 /* -------------------------------------------------------------------------- */ 1656 /* 1657 PCApplyTranspose_BDDC - Applies the transpose of the BDDC operator to a vector. 1658 1659 Input Parameters: 1660 + pc - the preconditioner context 1661 - r - input vector (global) 1662 1663 Output Parameter: 1664 . z - output vector (global) 1665 1666 Application Interface Routine: PCApplyTranspose() 1667 */ 1668 #undef __FUNCT__ 1669 #define __FUNCT__ "PCApplyTranspose_BDDC" 1670 PetscErrorCode PCApplyTranspose_BDDC(PC pc,Vec r,Vec z) 1671 { 1672 PC_IS *pcis = (PC_IS*)(pc->data); 1673 PC_BDDC *pcbddc = (PC_BDDC*)(pc->data); 1674 PetscInt n_B = pcis->n_B, n_D = pcis->n - n_B; 1675 PetscErrorCode ierr; 1676 const PetscScalar one = 1.0; 1677 const PetscScalar m_one = -1.0; 1678 const PetscScalar zero = 0.0; 1679 1680 PetscFunctionBegin; 1681 if (pcbddc->benign_saddle_point) { /* get p0 from r */ 1682 ierr = PCBDDCBenignGetOrSetP0(pc,r,PETSC_TRUE);CHKERRQ(ierr); 1683 } 1684 if (!pcbddc->use_exact_dirichlet_trick) { 1685 ierr = VecCopy(r,z);CHKERRQ(ierr); 1686 /* First Dirichlet solve */ 1687 ierr = VecScatterBegin(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1688 ierr = VecScatterEnd(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1689 /* 1690 Assembling right hand side for BDDC operator 1691 - pcis->vec1_D for the Dirichlet part (if needed, i.e. pcbddc->switch_static == PETSC_TRUE) 1692 - pcis->vec1_B the interface part of the global vector z 1693 */ 1694 if (n_D) { 1695 ierr = KSPSolveTranspose(pcbddc->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr); 1696 ierr = VecScale(pcis->vec2_D,m_one);CHKERRQ(ierr); 1697 if (pcbddc->switch_static) { ierr = MatMultTransposeAdd(pcis->A_II,pcis->vec2_D,pcis->vec1_D,pcis->vec1_D);CHKERRQ(ierr); } 1698 ierr = MatMultTranspose(pcis->A_IB,pcis->vec2_D,pcis->vec1_B);CHKERRQ(ierr); 1699 } else { 1700 ierr = VecSet(pcis->vec1_B,zero);CHKERRQ(ierr); 1701 } 1702 ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1703 ierr = VecScatterEnd(pcis->global_to_B,pcis->vec1_B,z,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1704 ierr = PCBDDCScalingRestriction(pc,z,pcis->vec1_B);CHKERRQ(ierr); 1705 } else { 1706 if (pcbddc->switch_static) { 1707 ierr = VecSet(pcis->vec1_D,zero);CHKERRQ(ierr); 1708 } 1709 ierr = PCBDDCScalingRestriction(pc,r,pcis->vec1_B);CHKERRQ(ierr); 1710 } 1711 1712 /* Apply interface preconditioner 1713 input/output vecs: pcis->vec1_B and pcis->vec1_D */ 1714 ierr = PCBDDCApplyInterfacePreconditioner(pc,PETSC_TRUE);CHKERRQ(ierr); 1715 1716 /* Apply transpose of partition of unity operator */ 1717 ierr = PCBDDCScalingExtension(pc,pcis->vec1_B,z);CHKERRQ(ierr); 1718 1719 /* Second Dirichlet solve and assembling of output */ 1720 ierr = VecScatterBegin(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1721 ierr = VecScatterEnd(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1722 if (n_B) { 1723 ierr = MatMultTranspose(pcis->A_BI,pcis->vec1_B,pcis->vec3_D);CHKERRQ(ierr); 1724 if (pcbddc->switch_static) { ierr = MatMultTransposeAdd(pcis->A_II,pcis->vec1_D,pcis->vec3_D,pcis->vec3_D);CHKERRQ(ierr); } 1725 } else if (pcbddc->switch_static) { 1726 ierr = MatMultTranspose(pcis->A_II,pcis->vec1_D,pcis->vec3_D);CHKERRQ(ierr); 1727 } 1728 ierr = KSPSolveTranspose(pcbddc->ksp_D,pcis->vec3_D,pcis->vec4_D);CHKERRQ(ierr); 1729 if (!pcbddc->use_exact_dirichlet_trick) { 1730 if (pcbddc->switch_static) { 1731 ierr = VecAXPBYPCZ(pcis->vec2_D,m_one,one,m_one,pcis->vec4_D,pcis->vec1_D);CHKERRQ(ierr); 1732 } else { 1733 ierr = VecAXPBY(pcis->vec2_D,m_one,m_one,pcis->vec4_D);CHKERRQ(ierr); 1734 } 1735 ierr = VecScatterBegin(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1736 ierr = VecScatterEnd(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1737 } else { 1738 if (pcbddc->switch_static) { 1739 ierr = VecAXPBY(pcis->vec4_D,one,m_one,pcis->vec1_D);CHKERRQ(ierr); 1740 } else { 1741 ierr = VecScale(pcis->vec4_D,m_one);CHKERRQ(ierr); 1742 } 1743 ierr = VecScatterBegin(pcis->global_to_D,pcis->vec4_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1744 ierr = VecScatterEnd(pcis->global_to_D,pcis->vec4_D,z,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1745 } 1746 if (pcbddc->benign_saddle_point) { /* set p0 (computed in PCBDDCApplyInterface) */ 1747 ierr = PCBDDCBenignGetOrSetP0(pc,z,PETSC_FALSE);CHKERRQ(ierr); 1748 } 1749 PetscFunctionReturn(0); 1750 } 1751 /* -------------------------------------------------------------------------- */ 1752 1753 #undef __FUNCT__ 1754 #define __FUNCT__ "PCDestroy_BDDC" 1755 PetscErrorCode PCDestroy_BDDC(PC pc) 1756 { 1757 PC_BDDC *pcbddc = (PC_BDDC*)pc->data; 1758 PetscErrorCode ierr; 1759 1760 PetscFunctionBegin; 1761 /* free data created by PCIS */ 1762 ierr = PCISDestroy(pc);CHKERRQ(ierr); 1763 /* free BDDC custom data */ 1764 ierr = PCBDDCResetCustomization(pc);CHKERRQ(ierr); 1765 /* destroy objects related to topography */ 1766 ierr = PCBDDCResetTopography(pc);CHKERRQ(ierr); 1767 /* free allocated graph structure */ 1768 ierr = PetscFree(pcbddc->mat_graph);CHKERRQ(ierr); 1769 /* free allocated sub schurs structure */ 1770 ierr = PetscFree(pcbddc->sub_schurs);CHKERRQ(ierr); 1771 /* destroy objects for scaling operator */ 1772 ierr = PCBDDCScalingDestroy(pc);CHKERRQ(ierr); 1773 ierr = PetscFree(pcbddc->deluxe_ctx);CHKERRQ(ierr); 1774 /* free solvers stuff */ 1775 ierr = PCBDDCResetSolvers(pc);CHKERRQ(ierr); 1776 /* free global vectors needed in presolve */ 1777 ierr = VecDestroy(&pcbddc->temp_solution);CHKERRQ(ierr); 1778 ierr = VecDestroy(&pcbddc->original_rhs);CHKERRQ(ierr); 1779 /* free stuff for change of basis hooks */ 1780 if (pcbddc->new_global_mat) { 1781 PCBDDCChange_ctx change_ctx; 1782 ierr = MatShellGetContext(pcbddc->new_global_mat,&change_ctx);CHKERRQ(ierr); 1783 ierr = MatDestroy(&change_ctx->original_mat);CHKERRQ(ierr); 1784 ierr = MatDestroy(&change_ctx->global_change);CHKERRQ(ierr); 1785 ierr = VecDestroyVecs(2,&change_ctx->work);CHKERRQ(ierr); 1786 ierr = PetscFree(change_ctx);CHKERRQ(ierr); 1787 } 1788 ierr = MatDestroy(&pcbddc->new_global_mat);CHKERRQ(ierr); 1789 /* remove functions */ 1790 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetChangeOfBasisMat_C",NULL);CHKERRQ(ierr); 1791 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",NULL);CHKERRQ(ierr); 1792 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",NULL);CHKERRQ(ierr); 1793 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",NULL);CHKERRQ(ierr); 1794 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",NULL);CHKERRQ(ierr); 1795 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",NULL);CHKERRQ(ierr); 1796 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",NULL);CHKERRQ(ierr); 1797 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",NULL);CHKERRQ(ierr); 1798 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundariesLocal_C",NULL);CHKERRQ(ierr); 1799 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",NULL);CHKERRQ(ierr); 1800 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundariesLocal_C",NULL);CHKERRQ(ierr); 1801 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",NULL);CHKERRQ(ierr); 1802 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundariesLocal_C",NULL);CHKERRQ(ierr); 1803 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",NULL);CHKERRQ(ierr); 1804 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundariesLocal_C",NULL);CHKERRQ(ierr); 1805 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",NULL);CHKERRQ(ierr); 1806 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplittingLocal_C",NULL);CHKERRQ(ierr); 1807 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",NULL);CHKERRQ(ierr); 1808 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",NULL);CHKERRQ(ierr); 1809 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",NULL);CHKERRQ(ierr); 1810 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",NULL);CHKERRQ(ierr); 1811 /* Free the private data structure */ 1812 ierr = PetscFree(pc->data);CHKERRQ(ierr); 1813 PetscFunctionReturn(0); 1814 } 1815 /* -------------------------------------------------------------------------- */ 1816 1817 #undef __FUNCT__ 1818 #define __FUNCT__ "PCBDDCMatFETIDPGetRHS_BDDC" 1819 static PetscErrorCode PCBDDCMatFETIDPGetRHS_BDDC(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs) 1820 { 1821 FETIDPMat_ctx mat_ctx; 1822 Vec copy_standard_rhs; 1823 PC_IS* pcis; 1824 PC_BDDC* pcbddc; 1825 PetscErrorCode ierr; 1826 1827 PetscFunctionBegin; 1828 ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr); 1829 pcis = (PC_IS*)mat_ctx->pc->data; 1830 pcbddc = (PC_BDDC*)mat_ctx->pc->data; 1831 1832 /* 1833 change of basis for physical rhs if needed 1834 It also changes the rhs in case of dirichlet boundaries 1835 TODO: better management when FETIDP will have its own class 1836 */ 1837 ierr = VecDuplicate(standard_rhs,©_standard_rhs);CHKERRQ(ierr); 1838 ierr = VecCopy(standard_rhs,copy_standard_rhs);CHKERRQ(ierr); 1839 ierr = PCPreSolve_BDDC(mat_ctx->pc,NULL,copy_standard_rhs,NULL);CHKERRQ(ierr); 1840 /* store vectors for computation of fetidp final solution */ 1841 ierr = VecScatterBegin(pcis->global_to_D,copy_standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1842 ierr = VecScatterEnd(pcis->global_to_D,copy_standard_rhs,mat_ctx->temp_solution_D,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1843 /* scale rhs since it should be unassembled */ 1844 /* TODO use counter scaling? (also below) */ 1845 ierr = VecScatterBegin(pcis->global_to_B,copy_standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1846 ierr = VecScatterEnd(pcis->global_to_B,copy_standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1847 /* Apply partition of unity */ 1848 ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr); 1849 /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,copy_standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */ 1850 if (!pcbddc->switch_static) { 1851 /* compute partially subassembled Schur complement right-hand side */ 1852 ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr); 1853 ierr = MatMult(pcis->A_BI,pcis->vec1_D,pcis->vec1_B);CHKERRQ(ierr); 1854 ierr = VecAXPY(mat_ctx->temp_solution_B,-1.0,pcis->vec1_B);CHKERRQ(ierr); 1855 ierr = VecSet(copy_standard_rhs,0.0);CHKERRQ(ierr); 1856 ierr = VecScatterBegin(pcis->global_to_B,mat_ctx->temp_solution_B,copy_standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1857 ierr = VecScatterEnd(pcis->global_to_B,mat_ctx->temp_solution_B,copy_standard_rhs,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1858 /* ierr = PCBDDCScalingRestriction(mat_ctx->pc,copy_standard_rhs,mat_ctx->temp_solution_B);CHKERRQ(ierr); */ 1859 ierr = VecScatterBegin(pcis->global_to_B,copy_standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1860 ierr = VecScatterEnd(pcis->global_to_B,copy_standard_rhs,mat_ctx->temp_solution_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1861 ierr = VecPointwiseMult(mat_ctx->temp_solution_B,pcis->D,mat_ctx->temp_solution_B);CHKERRQ(ierr); 1862 } 1863 ierr = VecDestroy(©_standard_rhs);CHKERRQ(ierr); 1864 /* BDDC rhs */ 1865 ierr = VecCopy(mat_ctx->temp_solution_B,pcis->vec1_B);CHKERRQ(ierr); 1866 if (pcbddc->switch_static) { 1867 ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr); 1868 } 1869 /* apply BDDC */ 1870 ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_FALSE);CHKERRQ(ierr); 1871 /* Application of B_delta and assembling of rhs for fetidp fluxes */ 1872 ierr = VecSet(fetidp_flux_rhs,0.0);CHKERRQ(ierr); 1873 ierr = MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local);CHKERRQ(ierr); 1874 ierr = VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1875 ierr = VecScatterEnd(mat_ctx->l2g_lambda,mat_ctx->lambda_local,fetidp_flux_rhs,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1876 PetscFunctionReturn(0); 1877 } 1878 1879 #undef __FUNCT__ 1880 #define __FUNCT__ "PCBDDCMatFETIDPGetRHS" 1881 /*@ 1882 PCBDDCMatFETIDPGetRHS - Compute the right-hand side for FETI-DP linear system using the physical right-hand side 1883 1884 Collective 1885 1886 Input Parameters: 1887 + fetidp_mat - the FETI-DP matrix object obtained by a call to PCBDDCCreateFETIDPOperators 1888 - standard_rhs - the right-hand side of the original linear system 1889 1890 Output Parameters: 1891 . fetidp_flux_rhs - the right-hand side for the FETI-DP linear system 1892 1893 Level: developer 1894 1895 Notes: 1896 1897 .seealso: PCBDDC, PCBDDCCreateFETIDPOperators, PCBDDCMatFETIDPGetSolution 1898 @*/ 1899 PetscErrorCode PCBDDCMatFETIDPGetRHS(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs) 1900 { 1901 FETIDPMat_ctx mat_ctx; 1902 PetscErrorCode ierr; 1903 1904 PetscFunctionBegin; 1905 ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr); 1906 ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetRHS_C",(Mat,Vec,Vec),(fetidp_mat,standard_rhs,fetidp_flux_rhs));CHKERRQ(ierr); 1907 PetscFunctionReturn(0); 1908 } 1909 /* -------------------------------------------------------------------------- */ 1910 1911 #undef __FUNCT__ 1912 #define __FUNCT__ "PCBDDCMatFETIDPGetSolution_BDDC" 1913 static PetscErrorCode PCBDDCMatFETIDPGetSolution_BDDC(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol) 1914 { 1915 FETIDPMat_ctx mat_ctx; 1916 PC_IS* pcis; 1917 PC_BDDC* pcbddc; 1918 PetscErrorCode ierr; 1919 1920 PetscFunctionBegin; 1921 ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr); 1922 pcis = (PC_IS*)mat_ctx->pc->data; 1923 pcbddc = (PC_BDDC*)mat_ctx->pc->data; 1924 1925 /* apply B_delta^T */ 1926 ierr = VecScatterBegin(mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1927 ierr = VecScatterEnd (mat_ctx->l2g_lambda,fetidp_flux_sol,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1928 ierr = MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B);CHKERRQ(ierr); 1929 /* compute rhs for BDDC application */ 1930 ierr = VecAYPX(pcis->vec1_B,-1.0,mat_ctx->temp_solution_B);CHKERRQ(ierr); 1931 if (pcbddc->switch_static) { 1932 ierr = VecCopy(mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr); 1933 } 1934 /* apply BDDC */ 1935 ierr = PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,PETSC_FALSE);CHKERRQ(ierr); 1936 /* put values into standard global vector */ 1937 ierr = VecScatterBegin(pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1938 ierr = VecScatterEnd (pcis->global_to_B,pcis->vec1_B,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1939 if (!pcbddc->switch_static) { 1940 /* compute values into the interior if solved for the partially subassembled Schur complement */ 1941 ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec1_D);CHKERRQ(ierr); 1942 ierr = VecAXPY(mat_ctx->temp_solution_D,-1.0,pcis->vec1_D);CHKERRQ(ierr); 1943 ierr = KSPSolve(pcbddc->ksp_D,mat_ctx->temp_solution_D,pcis->vec1_D);CHKERRQ(ierr); 1944 } 1945 ierr = VecScatterBegin(pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1946 ierr = VecScatterEnd (pcis->global_to_D,pcis->vec1_D,standard_sol,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1947 /* final change of basis if needed 1948 Is also sums the dirichlet part removed during RHS assembling */ 1949 ierr = PCPostSolve_BDDC(mat_ctx->pc,NULL,NULL,standard_sol);CHKERRQ(ierr); 1950 PetscFunctionReturn(0); 1951 } 1952 1953 #undef __FUNCT__ 1954 #define __FUNCT__ "PCBDDCMatFETIDPGetSolution" 1955 /*@ 1956 PCBDDCMatFETIDPGetSolution - Compute the physical solution using the solution of the FETI-DP linear system 1957 1958 Collective 1959 1960 Input Parameters: 1961 + fetidp_mat - the FETI-DP matrix obtained by a call to PCBDDCCreateFETIDPOperators 1962 - fetidp_flux_sol - the solution of the FETI-DP linear system 1963 1964 Output Parameters: 1965 . standard_sol - the solution defined on the physical domain 1966 1967 Level: developer 1968 1969 Notes: 1970 1971 .seealso: PCBDDC, PCBDDCCreateFETIDPOperators, PCBDDCMatFETIDPGetRHS 1972 @*/ 1973 PetscErrorCode PCBDDCMatFETIDPGetSolution(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol) 1974 { 1975 FETIDPMat_ctx mat_ctx; 1976 PetscErrorCode ierr; 1977 1978 PetscFunctionBegin; 1979 ierr = MatShellGetContext(fetidp_mat,&mat_ctx);CHKERRQ(ierr); 1980 ierr = PetscTryMethod(mat_ctx->pc,"PCBDDCMatFETIDPGetSolution_C",(Mat,Vec,Vec),(fetidp_mat,fetidp_flux_sol,standard_sol));CHKERRQ(ierr); 1981 PetscFunctionReturn(0); 1982 } 1983 /* -------------------------------------------------------------------------- */ 1984 1985 extern PetscErrorCode FETIDPMatMult(Mat,Vec,Vec); 1986 extern PetscErrorCode FETIDPMatMultTranspose(Mat,Vec,Vec); 1987 extern PetscErrorCode PCBDDCDestroyFETIDPMat(Mat); 1988 extern PetscErrorCode FETIDPPCApply(PC,Vec,Vec); 1989 extern PetscErrorCode FETIDPPCApplyTranspose(PC,Vec,Vec); 1990 extern PetscErrorCode PCBDDCDestroyFETIDPPC(PC); 1991 1992 #undef __FUNCT__ 1993 #define __FUNCT__ "PCBDDCCreateFETIDPOperators_BDDC" 1994 static PetscErrorCode PCBDDCCreateFETIDPOperators_BDDC(PC pc, Mat *fetidp_mat, PC *fetidp_pc) 1995 { 1996 1997 FETIDPMat_ctx fetidpmat_ctx; 1998 Mat newmat; 1999 FETIDPPC_ctx fetidppc_ctx; 2000 PC newpc; 2001 MPI_Comm comm; 2002 PetscErrorCode ierr; 2003 2004 PetscFunctionBegin; 2005 ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 2006 /* FETIDP linear matrix */ 2007 ierr = PCBDDCCreateFETIDPMatContext(pc,&fetidpmat_ctx);CHKERRQ(ierr); 2008 ierr = PCBDDCSetupFETIDPMatContext(fetidpmat_ctx);CHKERRQ(ierr); 2009 ierr = MatCreateShell(comm,PETSC_DECIDE,PETSC_DECIDE,fetidpmat_ctx->n_lambda,fetidpmat_ctx->n_lambda,fetidpmat_ctx,&newmat);CHKERRQ(ierr); 2010 ierr = MatShellSetOperation(newmat,MATOP_MULT,(void (*)(void))FETIDPMatMult);CHKERRQ(ierr); 2011 ierr = MatShellSetOperation(newmat,MATOP_MULT_TRANSPOSE,(void (*)(void))FETIDPMatMultTranspose);CHKERRQ(ierr); 2012 ierr = MatShellSetOperation(newmat,MATOP_DESTROY,(void (*)(void))PCBDDCDestroyFETIDPMat);CHKERRQ(ierr); 2013 ierr = MatSetUp(newmat);CHKERRQ(ierr); 2014 /* FETIDP preconditioner */ 2015 ierr = PCBDDCCreateFETIDPPCContext(pc,&fetidppc_ctx);CHKERRQ(ierr); 2016 ierr = PCBDDCSetupFETIDPPCContext(newmat,fetidppc_ctx);CHKERRQ(ierr); 2017 ierr = PCCreate(comm,&newpc);CHKERRQ(ierr); 2018 ierr = PCSetType(newpc,PCSHELL);CHKERRQ(ierr); 2019 ierr = PCShellSetContext(newpc,fetidppc_ctx);CHKERRQ(ierr); 2020 ierr = PCShellSetApply(newpc,FETIDPPCApply);CHKERRQ(ierr); 2021 ierr = PCShellSetApplyTranspose(newpc,FETIDPPCApplyTranspose);CHKERRQ(ierr); 2022 ierr = PCShellSetDestroy(newpc,PCBDDCDestroyFETIDPPC);CHKERRQ(ierr); 2023 ierr = PCSetOperators(newpc,newmat,newmat);CHKERRQ(ierr); 2024 ierr = PCSetUp(newpc);CHKERRQ(ierr); 2025 /* return pointers for objects created */ 2026 *fetidp_mat=newmat; 2027 *fetidp_pc=newpc; 2028 PetscFunctionReturn(0); 2029 } 2030 2031 #undef __FUNCT__ 2032 #define __FUNCT__ "PCBDDCCreateFETIDPOperators" 2033 /*@ 2034 PCBDDCCreateFETIDPOperators - Create FETI-DP operators 2035 2036 Collective 2037 2038 Input Parameters: 2039 . pc - the BDDC preconditioning context (setup should have been called before) 2040 2041 Output Parameters: 2042 + fetidp_mat - shell FETI-DP matrix object 2043 - fetidp_pc - shell Dirichlet preconditioner for FETI-DP matrix 2044 2045 Options Database Keys: 2046 . -fetidp_fullyredundant <false> - use or not a fully redundant set of Lagrange multipliers 2047 2048 Level: developer 2049 2050 Notes: 2051 Currently the only operations provided for FETI-DP matrix are MatMult and MatMultTranspose 2052 2053 .seealso: PCBDDC, PCBDDCMatFETIDPGetRHS, PCBDDCMatFETIDPGetSolution 2054 @*/ 2055 PetscErrorCode PCBDDCCreateFETIDPOperators(PC pc, Mat *fetidp_mat, PC *fetidp_pc) 2056 { 2057 PetscErrorCode ierr; 2058 2059 PetscFunctionBegin; 2060 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 2061 if (pc->setupcalled) { 2062 ierr = PetscUseMethod(pc,"PCBDDCCreateFETIDPOperators_C",(PC,Mat*,PC*),(pc,fetidp_mat,fetidp_pc));CHKERRQ(ierr); 2063 } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"You must call PCSetup_BDDC() first \n"); 2064 PetscFunctionReturn(0); 2065 } 2066 /* -------------------------------------------------------------------------- */ 2067 /*MC 2068 PCBDDC - Balancing Domain Decomposition by Constraints. 2069 2070 An implementation of the BDDC preconditioner based on 2071 2072 .vb 2073 [1] C. R. Dohrmann. "An approximate BDDC preconditioner", Numerical Linear Algebra with Applications Volume 14, Issue 2, pages 149-168, March 2007 2074 [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 2075 [3] J. Mandel, B. Sousedik, C. R. Dohrmann. "Multispace and Multilevel BDDC", http://arxiv.org/abs/0712.3977 2076 [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 2077 .ve 2078 2079 The matrix to be preconditioned (Pmat) must be of type MATIS. 2080 2081 Currently works with MATIS matrices with local matrices of type MATSEQAIJ, MATSEQBAIJ or MATSEQSBAIJ, either with real or complex numbers. 2082 2083 It also works with unsymmetric and indefinite problems. 2084 2085 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. 2086 2087 Approximate local solvers are automatically adapted for singular linear problems (see [1]) if the user has provided the nullspace using PCBDDCSetNullSpace() 2088 2089 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() 2090 Additional information on dofs can be provided by using PCBDDCSetDofsSplitting(), PCBDDCSetDirichletBoundaries(), PCBDDCSetNeumannBoundaries(), and PCBDDCSetPrimalVerticesLocalIS() 2091 2092 Constraints can be customized by attaching a MatNullSpace object to the MATIS matrix via MatSetNearNullSpace(). Non-singular modes are retained via SVD. 2093 2094 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. 2095 User defined change of basis can be passed to PCBDDC by using PCBDDCSetChangeOfBasisMat() 2096 2097 The PETSc implementation also supports multilevel BDDC [3]. Coarse grids are partitioned using a MatPartitioning object. 2098 2099 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. 2100 2101 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. 2102 Deluxe scaling is not supported yet for FETI-DP. 2103 2104 Options Database Keys (some of them, run with -h for a complete list): 2105 2106 . -pc_bddc_use_vertices <true> - use or not vertices in primal space 2107 . -pc_bddc_use_edges <true> - use or not edges in primal space 2108 . -pc_bddc_use_faces <false> - use or not faces in primal space 2109 . -pc_bddc_symmetric <true> - symmetric computation of primal basis functions. Specify false for unsymmetric problems 2110 . -pc_bddc_use_change_of_basis <false> - use change of basis approach (on edges only) 2111 . -pc_bddc_use_change_on_faces <false> - use change of basis approach on faces if change of basis has been requested 2112 . -pc_bddc_switch_static <false> - switches from M_2 (default) to M_3 operator (see reference article [1]) 2113 . -pc_bddc_levels <0> - maximum number of levels for multilevel 2114 . -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) 2115 . -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) 2116 . -pc_bddc_use_deluxe_scaling <false> - use deluxe scaling 2117 . -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) 2118 . -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) 2119 - -pc_bddc_check_level <0> - set verbosity level of debugging output 2120 2121 Options for Dirichlet, Neumann or coarse solver can be set with 2122 .vb 2123 -pc_bddc_dirichlet_ 2124 -pc_bddc_neumann_ 2125 -pc_bddc_coarse_ 2126 .ve 2127 e.g -pc_bddc_dirichlet_ksp_type richardson -pc_bddc_dirichlet_pc_type gamg. PCBDDC uses by default KPSPREONLY and PCLU. 2128 2129 When using a multilevel approach, solvers' options at the N-th level (N > 1) can be specified as 2130 .vb 2131 -pc_bddc_dirichlet_lN_ 2132 -pc_bddc_neumann_lN_ 2133 -pc_bddc_coarse_lN_ 2134 .ve 2135 Note that level number ranges from the finest (0) to the coarsest (N). 2136 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. 2137 .vb 2138 -pc_bddc_coarse_pc_bddc_adaptive_threshold 5 -pc_bddc_coarse_l1_pc_bddc_redistribute 3 2139 .ve 2140 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 2141 2142 Level: intermediate 2143 2144 Developer notes: 2145 2146 Contributed by Stefano Zampini 2147 2148 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, MATIS 2149 M*/ 2150 2151 #undef __FUNCT__ 2152 #define __FUNCT__ "PCCreate_BDDC" 2153 PETSC_EXTERN PetscErrorCode PCCreate_BDDC(PC pc) 2154 { 2155 PetscErrorCode ierr; 2156 PC_BDDC *pcbddc; 2157 2158 PetscFunctionBegin; 2159 /* Creates the private data structure for this preconditioner and attach it to the PC object. */ 2160 ierr = PetscNewLog(pc,&pcbddc);CHKERRQ(ierr); 2161 pc->data = (void*)pcbddc; 2162 2163 /* create PCIS data structure */ 2164 ierr = PCISCreate(pc);CHKERRQ(ierr); 2165 2166 /* BDDC customization */ 2167 pcbddc->use_local_adj = PETSC_TRUE; 2168 pcbddc->use_vertices = PETSC_TRUE; 2169 pcbddc->use_edges = PETSC_TRUE; 2170 pcbddc->use_faces = PETSC_FALSE; 2171 pcbddc->use_change_of_basis = PETSC_FALSE; 2172 pcbddc->use_change_on_faces = PETSC_FALSE; 2173 pcbddc->switch_static = PETSC_FALSE; 2174 pcbddc->use_nnsp_true = PETSC_FALSE; 2175 pcbddc->use_qr_single = PETSC_FALSE; 2176 pcbddc->symmetric_primal = PETSC_TRUE; 2177 pcbddc->benign_saddle_point = PETSC_FALSE; 2178 pcbddc->dbg_flag = 0; 2179 /* private */ 2180 pcbddc->local_primal_size = 0; 2181 pcbddc->local_primal_size_cc = 0; 2182 pcbddc->local_primal_ref_node = 0; 2183 pcbddc->local_primal_ref_mult = 0; 2184 pcbddc->n_vertices = 0; 2185 pcbddc->primal_indices_local_idxs = 0; 2186 pcbddc->recompute_topography = PETSC_FALSE; 2187 pcbddc->coarse_size = -1; 2188 pcbddc->new_primal_space = PETSC_FALSE; 2189 pcbddc->new_primal_space_local = PETSC_FALSE; 2190 pcbddc->global_primal_indices = 0; 2191 pcbddc->onearnullspace = 0; 2192 pcbddc->onearnullvecs_state = 0; 2193 pcbddc->user_primal_vertices = 0; 2194 pcbddc->NullSpace = 0; 2195 pcbddc->temp_solution = 0; 2196 pcbddc->original_rhs = 0; 2197 pcbddc->local_mat = 0; 2198 pcbddc->ChangeOfBasisMatrix = 0; 2199 pcbddc->user_ChangeOfBasisMatrix = 0; 2200 pcbddc->new_global_mat = 0; 2201 pcbddc->coarse_vec = 0; 2202 pcbddc->coarse_ksp = 0; 2203 pcbddc->coarse_phi_B = 0; 2204 pcbddc->coarse_phi_D = 0; 2205 pcbddc->coarse_psi_B = 0; 2206 pcbddc->coarse_psi_D = 0; 2207 pcbddc->vec1_P = 0; 2208 pcbddc->vec1_R = 0; 2209 pcbddc->vec2_R = 0; 2210 pcbddc->local_auxmat1 = 0; 2211 pcbddc->local_auxmat2 = 0; 2212 pcbddc->R_to_B = 0; 2213 pcbddc->R_to_D = 0; 2214 pcbddc->ksp_D = 0; 2215 pcbddc->ksp_R = 0; 2216 pcbddc->NeumannBoundaries = 0; 2217 pcbddc->NeumannBoundariesLocal = 0; 2218 pcbddc->DirichletBoundaries = 0; 2219 pcbddc->DirichletBoundariesLocal = 0; 2220 pcbddc->user_provided_isfordofs = PETSC_FALSE; 2221 pcbddc->n_ISForDofs = 0; 2222 pcbddc->n_ISForDofsLocal = 0; 2223 pcbddc->ISForDofs = 0; 2224 pcbddc->ISForDofsLocal = 0; 2225 pcbddc->ConstraintMatrix = 0; 2226 pcbddc->use_exact_dirichlet_trick = PETSC_TRUE; 2227 pcbddc->coarse_loc_to_glob = 0; 2228 pcbddc->coarsening_ratio = 8; 2229 pcbddc->coarse_adj_red = 0; 2230 pcbddc->current_level = 0; 2231 pcbddc->max_levels = 0; 2232 pcbddc->use_coarse_estimates = PETSC_FALSE; 2233 pcbddc->redistribute_coarse = 0; 2234 pcbddc->coarse_subassembling = 0; 2235 pcbddc->coarse_subassembling_init = 0; 2236 2237 /* benign subspace trick */ 2238 pcbddc->B0_ncol = 0; 2239 pcbddc->B0_cols = NULL; 2240 pcbddc->B0_vals = NULL; 2241 pcbddc->benign_change = 0; 2242 pcbddc->benign_vec = 0; 2243 pcbddc->benign_original_mat = 0; 2244 pcbddc->benign_sf = 0; 2245 pcbddc->benign_p0_lidx = -1; 2246 pcbddc->benign_p0_gidx = -1; 2247 pcbddc->benign_null = PETSC_FALSE; 2248 2249 /* create local graph structure */ 2250 ierr = PCBDDCGraphCreate(&pcbddc->mat_graph);CHKERRQ(ierr); 2251 2252 /* scaling */ 2253 pcbddc->work_scaling = 0; 2254 pcbddc->use_deluxe_scaling = PETSC_FALSE; 2255 pcbddc->faster_deluxe = PETSC_FALSE; 2256 2257 /* create sub schurs structure */ 2258 ierr = PCBDDCSubSchursCreate(&pcbddc->sub_schurs);CHKERRQ(ierr); 2259 pcbddc->sub_schurs_rebuild = PETSC_FALSE; 2260 pcbddc->sub_schurs_layers = -1; 2261 pcbddc->sub_schurs_use_useradj = PETSC_FALSE; 2262 2263 pcbddc->computed_rowadj = PETSC_FALSE; 2264 2265 /* adaptivity */ 2266 pcbddc->adaptive_threshold = 0.0; 2267 pcbddc->adaptive_nmax = 0; 2268 pcbddc->adaptive_nmin = 0; 2269 2270 /* function pointers */ 2271 pc->ops->apply = PCApply_BDDC; 2272 pc->ops->applytranspose = PCApplyTranspose_BDDC; 2273 pc->ops->setup = PCSetUp_BDDC; 2274 pc->ops->destroy = PCDestroy_BDDC; 2275 pc->ops->setfromoptions = PCSetFromOptions_BDDC; 2276 pc->ops->view = 0; 2277 pc->ops->applyrichardson = 0; 2278 pc->ops->applysymmetricleft = 0; 2279 pc->ops->applysymmetricright = 0; 2280 pc->ops->presolve = PCPreSolve_BDDC; 2281 pc->ops->postsolve = PCPostSolve_BDDC; 2282 2283 /* composing function */ 2284 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetChangeOfBasisMat_C",PCBDDCSetChangeOfBasisMat_BDDC);CHKERRQ(ierr); 2285 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetPrimalVerticesLocalIS_C",PCBDDCSetPrimalVerticesLocalIS_BDDC);CHKERRQ(ierr); 2286 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetCoarseningRatio_C",PCBDDCSetCoarseningRatio_BDDC);CHKERRQ(ierr); 2287 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevel_C",PCBDDCSetLevel_BDDC);CHKERRQ(ierr); 2288 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetUseExactDirichlet_C",PCBDDCSetUseExactDirichlet_BDDC);CHKERRQ(ierr); 2289 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLevels_C",PCBDDCSetLevels_BDDC);CHKERRQ(ierr); 2290 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNullSpace_C",PCBDDCSetNullSpace_BDDC);CHKERRQ(ierr); 2291 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundaries_C",PCBDDCSetDirichletBoundaries_BDDC);CHKERRQ(ierr); 2292 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDirichletBoundariesLocal_C",PCBDDCSetDirichletBoundariesLocal_BDDC);CHKERRQ(ierr); 2293 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundaries_C",PCBDDCSetNeumannBoundaries_BDDC);CHKERRQ(ierr); 2294 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetNeumannBoundariesLocal_C",PCBDDCSetNeumannBoundariesLocal_BDDC);CHKERRQ(ierr); 2295 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundaries_C",PCBDDCGetDirichletBoundaries_BDDC);CHKERRQ(ierr); 2296 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetDirichletBoundariesLocal_C",PCBDDCGetDirichletBoundariesLocal_BDDC);CHKERRQ(ierr); 2297 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundaries_C",PCBDDCGetNeumannBoundaries_BDDC);CHKERRQ(ierr); 2298 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCGetNeumannBoundariesLocal_C",PCBDDCGetNeumannBoundariesLocal_BDDC);CHKERRQ(ierr); 2299 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplitting_C",PCBDDCSetDofsSplitting_BDDC);CHKERRQ(ierr); 2300 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetDofsSplittingLocal_C",PCBDDCSetDofsSplittingLocal_BDDC);CHKERRQ(ierr); 2301 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCSetLocalAdjacencyGraph_C",PCBDDCSetLocalAdjacencyGraph_BDDC);CHKERRQ(ierr); 2302 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCCreateFETIDPOperators_C",PCBDDCCreateFETIDPOperators_BDDC);CHKERRQ(ierr); 2303 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetRHS_C",PCBDDCMatFETIDPGetRHS_BDDC);CHKERRQ(ierr); 2304 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCBDDCMatFETIDPGetSolution_C",PCBDDCMatFETIDPGetSolution_BDDC);CHKERRQ(ierr); 2305 PetscFunctionReturn(0); 2306 } 2307 2308