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