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