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 MATIS related operations contained in BDDC code 18 - Provide general case for subassembling 19 20 */ 21 22 #include <petsc/private/pcbddcimpl.h> /*I "petscpc.h" I*/ /* includes for fortran wrappers */ 23 #include <petsc/private/pcbddcprivateimpl.h> 24 #include <petscblaslapack.h> 25 26 static PetscBool PCBDDCPackageInitialized = PETSC_FALSE; 27 28 static PetscBool cited = PETSC_FALSE; 29 static const char citation[] = "@article{ZampiniPCBDDC,\n" 30 "author = {Stefano Zampini},\n" 31 "title = {{PCBDDC}: A Class of Robust Dual-Primal Methods in {PETS}c},\n" 32 "journal = {SIAM Journal on Scientific Computing},\n" 33 "volume = {38},\n" 34 "number = {5},\n" 35 "pages = {S282-S306},\n" 36 "year = {2016},\n" 37 "doi = {10.1137/15M1025785},\n" 38 "URL = {http://dx.doi.org/10.1137/15M1025785},\n" 39 "eprint = {http://dx.doi.org/10.1137/15M1025785}\n" 40 "}\n"; 41 42 PetscLogEvent PC_BDDC_Topology[PETSC_PCBDDC_MAXLEVELS]; 43 PetscLogEvent PC_BDDC_LocalSolvers[PETSC_PCBDDC_MAXLEVELS]; 44 PetscLogEvent PC_BDDC_LocalWork[PETSC_PCBDDC_MAXLEVELS]; 45 PetscLogEvent PC_BDDC_CorrectionSetUp[PETSC_PCBDDC_MAXLEVELS]; 46 PetscLogEvent PC_BDDC_ApproxSetUp[PETSC_PCBDDC_MAXLEVELS]; 47 PetscLogEvent PC_BDDC_ApproxApply[PETSC_PCBDDC_MAXLEVELS]; 48 PetscLogEvent PC_BDDC_CoarseSetUp[PETSC_PCBDDC_MAXLEVELS]; 49 PetscLogEvent PC_BDDC_CoarseSolver[PETSC_PCBDDC_MAXLEVELS]; 50 PetscLogEvent PC_BDDC_AdaptiveSetUp[PETSC_PCBDDC_MAXLEVELS]; 51 PetscLogEvent PC_BDDC_Scaling[PETSC_PCBDDC_MAXLEVELS]; 52 PetscLogEvent PC_BDDC_Schurs[PETSC_PCBDDC_MAXLEVELS]; 53 PetscLogEvent PC_BDDC_Solves[PETSC_PCBDDC_MAXLEVELS][3]; 54 55 const char *const PCBDDCInterfaceExtTypes[] = {"DIRICHLET", "LUMP", "PCBDDCInterfaceExtType", "PC_BDDC_INTERFACE_EXT_", NULL}; 56 57 PetscErrorCode PCApply_BDDC(PC, Vec, Vec); 58 59 PetscErrorCode PCSetFromOptions_BDDC(PC pc, PetscOptionItems *PetscOptionsObject) 60 { 61 PC_BDDC *pcbddc = (PC_BDDC *)pc->data; 62 PetscInt nt, i; 63 64 PetscFunctionBegin; 65 PetscOptionsHeadBegin(PetscOptionsObject, "BDDC options"); 66 /* Verbose debugging */ 67 PetscCall(PetscOptionsInt("-pc_bddc_check_level", "Verbose output for PCBDDC (intended for debug)", "none", pcbddc->dbg_flag, &pcbddc->dbg_flag, NULL)); 68 /* Approximate solvers */ 69 PetscCall(PetscOptionsEnum("-pc_bddc_interface_ext_type", "Use DIRICHLET or LUMP to extend interface corrections to interior", "PCBDDCSetInterfaceExtType", PCBDDCInterfaceExtTypes, (PetscEnum)pcbddc->interface_extension, (PetscEnum *)&pcbddc->interface_extension, NULL)); 70 if (pcbddc->interface_extension == PC_BDDC_INTERFACE_EXT_DIRICHLET) { 71 PetscCall(PetscOptionsBool("-pc_bddc_dirichlet_approximate", "Inform PCBDDC that we are using approximate Dirichlet solvers", "none", pcbddc->NullSpace_corr[0], &pcbddc->NullSpace_corr[0], NULL)); 72 PetscCall(PetscOptionsBool("-pc_bddc_dirichlet_approximate_scale", "Inform PCBDDC that we need to scale the Dirichlet solve", "none", pcbddc->NullSpace_corr[1], &pcbddc->NullSpace_corr[1], NULL)); 73 } else { 74 /* This flag is needed/implied by lumping */ 75 pcbddc->switch_static = PETSC_TRUE; 76 } 77 PetscCall(PetscOptionsBool("-pc_bddc_neumann_approximate", "Inform PCBDDC that we are using approximate Neumann solvers", "none", pcbddc->NullSpace_corr[2], &pcbddc->NullSpace_corr[2], NULL)); 78 PetscCall(PetscOptionsBool("-pc_bddc_neumann_approximate_scale", "Inform PCBDDC that we need to scale the Neumann solve", "none", pcbddc->NullSpace_corr[3], &pcbddc->NullSpace_corr[3], NULL)); 79 /* Primal space customization */ 80 PetscCall(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)); 81 PetscCall(PetscOptionsInt("-pc_bddc_graph_maxcount", "Maximum number of shared subdomains for a connected component", "none", pcbddc->graphmaxcount, &pcbddc->graphmaxcount, NULL)); 82 PetscCall(PetscOptionsBool("-pc_bddc_corner_selection", "Activates face-based corner selection", "none", pcbddc->corner_selection, &pcbddc->corner_selection, NULL)); 83 PetscCall(PetscOptionsBool("-pc_bddc_use_vertices", "Use or not corner dofs in coarse space", "none", pcbddc->use_vertices, &pcbddc->use_vertices, NULL)); 84 PetscCall(PetscOptionsBool("-pc_bddc_use_edges", "Use or not edge constraints in coarse space", "none", pcbddc->use_edges, &pcbddc->use_edges, NULL)); 85 PetscCall(PetscOptionsBool("-pc_bddc_use_faces", "Use or not face constraints in coarse space", "none", pcbddc->use_faces, &pcbddc->use_faces, NULL)); 86 PetscCall(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)); 87 PetscCall(PetscOptionsBool("-pc_bddc_use_nnsp", "Use near null space attached to the matrix to compute constraints", "none", pcbddc->use_nnsp, &pcbddc->use_nnsp, NULL)); 88 PetscCall(PetscOptionsBool("-pc_bddc_use_nnsp_true", "Use near null space attached to the matrix to compute constraints as is", "none", pcbddc->use_nnsp_true, &pcbddc->use_nnsp_true, NULL)); 89 PetscCall(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)); 90 /* Change of basis */ 91 PetscCall(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)); 92 PetscCall(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)); 93 if (!pcbddc->use_change_of_basis) pcbddc->use_change_on_faces = PETSC_FALSE; 94 /* Switch between M_2 (default) and M_3 preconditioners (as defined by C. Dohrmann in the ref. article) */ 95 PetscCall(PetscOptionsBool("-pc_bddc_switch_static", "Switch on static condensation ops around the interface preconditioner", "none", pcbddc->switch_static, &pcbddc->switch_static, NULL)); 96 PetscCall(PetscOptionsInt("-pc_bddc_coarse_eqs_per_proc", "Target number of equations per process for coarse problem redistribution (significant only at the coarsest level)", "none", pcbddc->coarse_eqs_per_proc, &pcbddc->coarse_eqs_per_proc, NULL)); 97 i = pcbddc->coarsening_ratio; 98 PetscCall(PetscOptionsInt("-pc_bddc_coarsening_ratio", "Set coarsening ratio used in multilevel coarsening", "PCBDDCSetCoarseningRatio", i, &i, NULL)); 99 PetscCall(PCBDDCSetCoarseningRatio(pc, i)); 100 i = pcbddc->max_levels; 101 PetscCall(PetscOptionsInt("-pc_bddc_levels", "Set maximum number of levels for multilevel", "PCBDDCSetLevels", i, &i, NULL)); 102 PetscCall(PCBDDCSetLevels(pc, i)); 103 PetscCall(PetscOptionsInt("-pc_bddc_coarse_eqs_limit", "Set maximum number of equations on coarsest grid to aim for", "none", pcbddc->coarse_eqs_limit, &pcbddc->coarse_eqs_limit, NULL)); 104 PetscCall(PetscOptionsBool("-pc_bddc_use_coarse_estimates", "Use estimated eigenvalues for coarse problem", "none", pcbddc->use_coarse_estimates, &pcbddc->use_coarse_estimates, NULL)); 105 PetscCall(PetscOptionsBool("-pc_bddc_use_deluxe_scaling", "Use deluxe scaling for BDDC", "none", pcbddc->use_deluxe_scaling, &pcbddc->use_deluxe_scaling, NULL)); 106 PetscCall(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)); 107 PetscCall(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)); 108 PetscCall(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)); 109 PetscCall(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)); 110 PetscCall(PetscOptionsBool("-pc_bddc_deluxe_zerorows", "Zero rows and columns of deluxe operators associated with primal dofs", "none", pcbddc->deluxe_zerorows, &pcbddc->deluxe_zerorows, NULL)); 111 PetscCall(PetscOptionsBool("-pc_bddc_deluxe_singlemat", "Collapse deluxe operators", "none", pcbddc->deluxe_singlemat, &pcbddc->deluxe_singlemat, NULL)); 112 PetscCall(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)); 113 nt = 2; 114 PetscCall(PetscOptionsRealArray("-pc_bddc_adaptive_threshold", "Thresholds to be used for adaptive selection of constraints", "none", pcbddc->adaptive_threshold, &nt, NULL)); 115 if (nt == 1) pcbddc->adaptive_threshold[1] = pcbddc->adaptive_threshold[0]; 116 PetscCall(PetscOptionsInt("-pc_bddc_adaptive_nmin", "Minimum number of constraints per connected components", "none", pcbddc->adaptive_nmin, &pcbddc->adaptive_nmin, NULL)); 117 PetscCall(PetscOptionsInt("-pc_bddc_adaptive_nmax", "Maximum number of constraints per connected components", "none", pcbddc->adaptive_nmax, &pcbddc->adaptive_nmax, NULL)); 118 PetscCall(PetscOptionsBool("-pc_bddc_symmetric", "Symmetric computation of primal basis functions", "none", pcbddc->symmetric_primal, &pcbddc->symmetric_primal, NULL)); 119 PetscCall(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)); 120 PetscCall(PetscOptionsBool("-pc_bddc_benign_trick", "Apply the benign subspace trick to saddle point problems with discontinuous pressures", "none", pcbddc->benign_saddle_point, &pcbddc->benign_saddle_point, NULL)); 121 PetscCall(PetscOptionsBool("-pc_bddc_benign_change", "Compute the pressure change of basis explicitly", "none", pcbddc->benign_change_explicit, &pcbddc->benign_change_explicit, NULL)); 122 PetscCall(PetscOptionsBool("-pc_bddc_benign_compute_correction", "Compute the benign correction during PreSolve", "none", pcbddc->benign_compute_correction, &pcbddc->benign_compute_correction, NULL)); 123 PetscCall(PetscOptionsBool("-pc_bddc_nonetflux", "Automatic computation of no-net-flux quadrature weights", "none", pcbddc->compute_nonetflux, &pcbddc->compute_nonetflux, NULL)); 124 PetscCall(PetscOptionsBool("-pc_bddc_detect_disconnected", "Detects disconnected subdomains", "none", pcbddc->detect_disconnected, &pcbddc->detect_disconnected, NULL)); 125 PetscCall(PetscOptionsBool("-pc_bddc_detect_disconnected_filter", "Filters out small entries in the local matrix when detecting disconnected subdomains", "none", pcbddc->detect_disconnected_filter, &pcbddc->detect_disconnected_filter, NULL)); 126 PetscCall(PetscOptionsBool("-pc_bddc_eliminate_dirichlet", "Whether or not we want to eliminate dirichlet dofs during presolve", "none", pcbddc->eliminate_dirdofs, &pcbddc->eliminate_dirdofs, NULL)); 127 PetscOptionsHeadEnd(); 128 PetscFunctionReturn(0); 129 } 130 131 static PetscErrorCode PCView_BDDC(PC pc, PetscViewer viewer) 132 { 133 PC_BDDC *pcbddc = (PC_BDDC *)pc->data; 134 PC_IS *pcis = (PC_IS *)pc->data; 135 PetscBool isascii; 136 PetscSubcomm subcomm; 137 PetscViewer subviewer; 138 139 PetscFunctionBegin; 140 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii)); 141 /* ASCII viewer */ 142 if (isascii) { 143 PetscMPIInt color, rank, size; 144 PetscInt64 loc[7], gsum[6], gmax[6], gmin[6], totbenign; 145 PetscScalar interface_size; 146 PetscReal ratio1 = 0., ratio2 = 0.; 147 Vec counter; 148 149 if (!pc->setupcalled) PetscCall(PetscViewerASCIIPrintf(viewer, " Partial information available: preconditioner has not been setup yet\n")); 150 PetscCall(PetscViewerASCIIPrintf(viewer, " Use verbose output: %" PetscInt_FMT "\n", pcbddc->dbg_flag)); 151 PetscCall(PetscViewerASCIIPrintf(viewer, " Use user-defined CSR: %d\n", !!pcbddc->mat_graph->nvtxs_csr)); 152 PetscCall(PetscViewerASCIIPrintf(viewer, " Use local mat graph: %d\n", pcbddc->use_local_adj && !pcbddc->mat_graph->nvtxs_csr)); 153 if (pcbddc->mat_graph->twodim) { 154 PetscCall(PetscViewerASCIIPrintf(viewer, " Connectivity graph topological dimension: 2\n")); 155 } else { 156 PetscCall(PetscViewerASCIIPrintf(viewer, " Connectivity graph topological dimension: 3\n")); 157 } 158 if (pcbddc->graphmaxcount != PETSC_MAX_INT) PetscCall(PetscViewerASCIIPrintf(viewer, " Graph max count: %" PetscInt_FMT "\n", pcbddc->graphmaxcount)); 159 PetscCall(PetscViewerASCIIPrintf(viewer, " Corner selection: %d (selected %d)\n", pcbddc->corner_selection, pcbddc->corner_selected)); 160 PetscCall(PetscViewerASCIIPrintf(viewer, " Use vertices: %d (vertex size %" PetscInt_FMT ")\n", pcbddc->use_vertices, pcbddc->vertex_size)); 161 PetscCall(PetscViewerASCIIPrintf(viewer, " Use edges: %d\n", pcbddc->use_edges)); 162 PetscCall(PetscViewerASCIIPrintf(viewer, " Use faces: %d\n", pcbddc->use_faces)); 163 PetscCall(PetscViewerASCIIPrintf(viewer, " Use true near null space: %d\n", pcbddc->use_nnsp_true)); 164 PetscCall(PetscViewerASCIIPrintf(viewer, " Use QR for single constraints on cc: %d\n", pcbddc->use_qr_single)); 165 PetscCall(PetscViewerASCIIPrintf(viewer, " Use change of basis on local edge nodes: %d\n", pcbddc->use_change_of_basis)); 166 PetscCall(PetscViewerASCIIPrintf(viewer, " Use change of basis on local face nodes: %d\n", pcbddc->use_change_on_faces)); 167 PetscCall(PetscViewerASCIIPrintf(viewer, " User defined change of basis matrix: %d\n", !!pcbddc->user_ChangeOfBasisMatrix)); 168 PetscCall(PetscViewerASCIIPrintf(viewer, " Has change of basis matrix: %d\n", !!pcbddc->ChangeOfBasisMatrix)); 169 PetscCall(PetscViewerASCIIPrintf(viewer, " Eliminate dirichlet boundary dofs: %d\n", pcbddc->eliminate_dirdofs)); 170 PetscCall(PetscViewerASCIIPrintf(viewer, " Switch on static condensation ops around the interface preconditioner: %d\n", pcbddc->switch_static)); 171 PetscCall(PetscViewerASCIIPrintf(viewer, " Use exact dirichlet trick: %d\n", pcbddc->use_exact_dirichlet_trick)); 172 PetscCall(PetscViewerASCIIPrintf(viewer, " Interface extension: %s\n", PCBDDCInterfaceExtTypes[pcbddc->interface_extension])); 173 PetscCall(PetscViewerASCIIPrintf(viewer, " Multilevel max levels: %" PetscInt_FMT "\n", pcbddc->max_levels)); 174 PetscCall(PetscViewerASCIIPrintf(viewer, " Multilevel coarsening ratio: %" PetscInt_FMT "\n", pcbddc->coarsening_ratio)); 175 PetscCall(PetscViewerASCIIPrintf(viewer, " Use estimated eigs for coarse problem: %d\n", pcbddc->use_coarse_estimates)); 176 PetscCall(PetscViewerASCIIPrintf(viewer, " Use deluxe scaling: %d\n", pcbddc->use_deluxe_scaling)); 177 PetscCall(PetscViewerASCIIPrintf(viewer, " Use deluxe zerorows: %d\n", pcbddc->deluxe_zerorows)); 178 PetscCall(PetscViewerASCIIPrintf(viewer, " Use deluxe singlemat: %d\n", pcbddc->deluxe_singlemat)); 179 PetscCall(PetscViewerASCIIPrintf(viewer, " Rebuild interface graph for Schur principal minors: %d\n", pcbddc->sub_schurs_rebuild)); 180 PetscCall(PetscViewerASCIIPrintf(viewer, " Number of dofs' layers for the computation of principal minors: %" PetscInt_FMT "\n", pcbddc->sub_schurs_layers)); 181 PetscCall(PetscViewerASCIIPrintf(viewer, " Use user CSR graph to compute successive layers: %d\n", pcbddc->sub_schurs_use_useradj)); 182 if (pcbddc->adaptive_threshold[1] != pcbddc->adaptive_threshold[0]) { 183 PetscCall(PetscViewerASCIIPrintf(viewer, " Adaptive constraint selection thresholds (active %d, userdefined %d): %g,%g\n", pcbddc->adaptive_selection, pcbddc->adaptive_userdefined, (double)pcbddc->adaptive_threshold[0], (double)pcbddc->adaptive_threshold[1])); 184 } else { 185 PetscCall(PetscViewerASCIIPrintf(viewer, " Adaptive constraint selection threshold (active %d, userdefined %d): %g\n", pcbddc->adaptive_selection, pcbddc->adaptive_userdefined, (double)pcbddc->adaptive_threshold[0])); 186 } 187 PetscCall(PetscViewerASCIIPrintf(viewer, " Min constraints / connected component: %" PetscInt_FMT "\n", pcbddc->adaptive_nmin)); 188 PetscCall(PetscViewerASCIIPrintf(viewer, " Max constraints / connected component: %" PetscInt_FMT "\n", pcbddc->adaptive_nmax)); 189 PetscCall(PetscViewerASCIIPrintf(viewer, " Invert exact Schur complement for adaptive selection: %d\n", pcbddc->sub_schurs_exact_schur)); 190 PetscCall(PetscViewerASCIIPrintf(viewer, " Symmetric computation of primal basis functions: %d\n", pcbddc->symmetric_primal)); 191 PetscCall(PetscViewerASCIIPrintf(viewer, " Num. Procs. to map coarse adjacency list: %" PetscInt_FMT "\n", pcbddc->coarse_adj_red)); 192 PetscCall(PetscViewerASCIIPrintf(viewer, " Coarse eqs per proc (significant at the coarsest level): %" PetscInt_FMT "\n", pcbddc->coarse_eqs_per_proc)); 193 PetscCall(PetscViewerASCIIPrintf(viewer, " Detect disconnected: %d (filter %d)\n", pcbddc->detect_disconnected, pcbddc->detect_disconnected_filter)); 194 PetscCall(PetscViewerASCIIPrintf(viewer, " Benign subspace trick: %d (change explicit %d)\n", pcbddc->benign_saddle_point, pcbddc->benign_change_explicit)); 195 PetscCall(PetscViewerASCIIPrintf(viewer, " Benign subspace trick is active: %d\n", pcbddc->benign_have_null)); 196 PetscCall(PetscViewerASCIIPrintf(viewer, " Algebraic computation of no-net-flux: %d\n", pcbddc->compute_nonetflux)); 197 if (!pc->setupcalled) PetscFunctionReturn(0); 198 199 /* compute interface size */ 200 PetscCall(VecSet(pcis->vec1_B, 1.0)); 201 PetscCall(MatCreateVecs(pc->pmat, &counter, NULL)); 202 PetscCall(VecSet(counter, 0.0)); 203 PetscCall(VecScatterBegin(pcis->global_to_B, pcis->vec1_B, counter, INSERT_VALUES, SCATTER_REVERSE)); 204 PetscCall(VecScatterEnd(pcis->global_to_B, pcis->vec1_B, counter, INSERT_VALUES, SCATTER_REVERSE)); 205 PetscCall(VecSum(counter, &interface_size)); 206 PetscCall(VecDestroy(&counter)); 207 208 /* compute some statistics on the domain decomposition */ 209 gsum[0] = 1; 210 gsum[1] = gsum[2] = gsum[3] = gsum[4] = gsum[5] = 0; 211 loc[0] = !!pcis->n; 212 loc[1] = pcis->n - pcis->n_B; 213 loc[2] = pcis->n_B; 214 loc[3] = pcbddc->local_primal_size; 215 loc[4] = pcis->n; 216 loc[5] = pcbddc->n_local_subs > 0 ? pcbddc->n_local_subs : (pcis->n ? 1 : 0); 217 loc[6] = pcbddc->benign_n; 218 PetscCallMPI(MPI_Reduce(loc, gsum, 6, MPIU_INT64, MPI_SUM, 0, PetscObjectComm((PetscObject)pc))); 219 if (!loc[0]) loc[1] = loc[2] = loc[3] = loc[4] = loc[5] = -1; 220 PetscCallMPI(MPI_Reduce(loc, gmax, 6, MPIU_INT64, MPI_MAX, 0, PetscObjectComm((PetscObject)pc))); 221 if (!loc[0]) loc[1] = loc[2] = loc[3] = loc[4] = loc[5] = PETSC_MAX_INT; 222 PetscCallMPI(MPI_Reduce(loc, gmin, 6, MPIU_INT64, MPI_MIN, 0, PetscObjectComm((PetscObject)pc))); 223 PetscCallMPI(MPI_Reduce(&loc[6], &totbenign, 1, MPIU_INT64, MPI_SUM, 0, PetscObjectComm((PetscObject)pc))); 224 if (pcbddc->coarse_size) { 225 ratio1 = pc->pmat->rmap->N / (1. * pcbddc->coarse_size); 226 ratio2 = PetscRealPart(interface_size) / pcbddc->coarse_size; 227 } 228 PetscCall(PetscViewerASCIIPrintf(viewer, "********************************** STATISTICS AT LEVEL %" PetscInt_FMT " **********************************\n", pcbddc->current_level)); 229 PetscCall(PetscViewerASCIIPrintf(viewer, " Global dofs sizes: all %" PetscInt_FMT " interface %" PetscInt_FMT " coarse %" PetscInt_FMT "\n", pc->pmat->rmap->N, (PetscInt)PetscRealPart(interface_size), pcbddc->coarse_size)); 230 PetscCall(PetscViewerASCIIPrintf(viewer, " Coarsening ratios: all/coarse %" PetscInt_FMT " interface/coarse %" PetscInt_FMT "\n", (PetscInt)ratio1, (PetscInt)ratio2)); 231 PetscCall(PetscViewerASCIIPrintf(viewer, " Active processes : %" PetscInt_FMT "\n", (PetscInt)gsum[0])); 232 PetscCall(PetscViewerASCIIPrintf(viewer, " Total subdomains : %" PetscInt_FMT "\n", (PetscInt)gsum[5])); 233 if (pcbddc->benign_have_null) PetscCall(PetscViewerASCIIPrintf(viewer, " Benign subs : %" PetscInt_FMT "\n", (PetscInt)totbenign)); 234 PetscCall(PetscViewerASCIIPrintf(viewer, " Dofs type :\tMIN\tMAX\tMEAN\n")); 235 PetscCall(PetscViewerASCIIPrintf(viewer, " Interior dofs :\t%" PetscInt_FMT "\t%" PetscInt_FMT "\t%" PetscInt_FMT "\n", (PetscInt)gmin[1], (PetscInt)gmax[1], (PetscInt)(gsum[1] / gsum[0]))); 236 PetscCall(PetscViewerASCIIPrintf(viewer, " Interface dofs :\t%" PetscInt_FMT "\t%" PetscInt_FMT "\t%" PetscInt_FMT "\n", (PetscInt)gmin[2], (PetscInt)gmax[2], (PetscInt)(gsum[2] / gsum[0]))); 237 PetscCall(PetscViewerASCIIPrintf(viewer, " Primal dofs :\t%" PetscInt_FMT "\t%" PetscInt_FMT "\t%" PetscInt_FMT "\n", (PetscInt)gmin[3], (PetscInt)gmax[3], (PetscInt)(gsum[3] / gsum[0]))); 238 PetscCall(PetscViewerASCIIPrintf(viewer, " Local dofs :\t%" PetscInt_FMT "\t%" PetscInt_FMT "\t%" PetscInt_FMT "\n", (PetscInt)gmin[4], (PetscInt)gmax[4], (PetscInt)(gsum[4] / gsum[0]))); 239 PetscCall(PetscViewerASCIIPrintf(viewer, " Local subs :\t%" PetscInt_FMT "\t%" PetscInt_FMT "\n", (PetscInt)gmin[5], (PetscInt)gmax[5])); 240 PetscCall(PetscViewerFlush(viewer)); 241 242 PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank)); 243 244 /* local solvers */ 245 PetscCall(PetscViewerGetSubViewer(viewer, PetscObjectComm((PetscObject)pcbddc->ksp_D), &subviewer)); 246 if (rank == 0) { 247 PetscCall(PetscViewerASCIIPrintf(subviewer, "--- Interior solver (rank 0)\n")); 248 PetscCall(PetscViewerASCIIPushTab(subviewer)); 249 PetscCall(KSPView(pcbddc->ksp_D, subviewer)); 250 PetscCall(PetscViewerASCIIPopTab(subviewer)); 251 PetscCall(PetscViewerASCIIPrintf(subviewer, "--- Correction solver (rank 0)\n")); 252 PetscCall(PetscViewerASCIIPushTab(subviewer)); 253 PetscCall(KSPView(pcbddc->ksp_R, subviewer)); 254 PetscCall(PetscViewerASCIIPopTab(subviewer)); 255 PetscCall(PetscViewerFlush(subviewer)); 256 } 257 PetscCall(PetscViewerRestoreSubViewer(viewer, PetscObjectComm((PetscObject)pcbddc->ksp_D), &subviewer)); 258 PetscCall(PetscViewerFlush(viewer)); 259 260 /* the coarse problem can be handled by a different communicator */ 261 if (pcbddc->coarse_ksp) color = 1; 262 else color = 0; 263 PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size)); 264 PetscCall(PetscSubcommCreate(PetscObjectComm((PetscObject)pc), &subcomm)); 265 PetscCall(PetscSubcommSetNumber(subcomm, PetscMin(size, 2))); 266 PetscCall(PetscSubcommSetTypeGeneral(subcomm, color, rank)); 267 PetscCall(PetscViewerGetSubViewer(viewer, PetscSubcommChild(subcomm), &subviewer)); 268 if (color == 1) { 269 PetscCall(PetscViewerASCIIPrintf(subviewer, "--- Coarse solver\n")); 270 PetscCall(PetscViewerASCIIPushTab(subviewer)); 271 PetscCall(KSPView(pcbddc->coarse_ksp, subviewer)); 272 PetscCall(PetscViewerASCIIPopTab(subviewer)); 273 PetscCall(PetscViewerFlush(subviewer)); 274 } 275 PetscCall(PetscViewerRestoreSubViewer(viewer, PetscSubcommChild(subcomm), &subviewer)); 276 PetscCall(PetscSubcommDestroy(&subcomm)); 277 PetscCall(PetscViewerFlush(viewer)); 278 } 279 PetscFunctionReturn(0); 280 } 281 282 static PetscErrorCode PCBDDCSetDiscreteGradient_BDDC(PC pc, Mat G, PetscInt order, PetscInt field, PetscBool global, PetscBool conforming) 283 { 284 PC_BDDC *pcbddc = (PC_BDDC *)pc->data; 285 286 PetscFunctionBegin; 287 PetscCall(PetscObjectReference((PetscObject)G)); 288 PetscCall(MatDestroy(&pcbddc->discretegradient)); 289 pcbddc->discretegradient = G; 290 pcbddc->nedorder = order > 0 ? order : -order; 291 pcbddc->nedfield = field; 292 pcbddc->nedglobal = global; 293 pcbddc->conforming = conforming; 294 PetscFunctionReturn(0); 295 } 296 297 /*@ 298 PCBDDCSetDiscreteGradient - Sets the discrete gradient 299 300 Collective on pc 301 302 Input Parameters: 303 + pc - the preconditioning context 304 . G - the discrete gradient matrix (in `MATAIJ` format) 305 . order - the order of the Nedelec space (1 for the lowest order) 306 . field - the field id of the Nedelec dofs (not used if the fields have not been specified) 307 . global - the type of global ordering for the rows of G 308 - conforming - whether the mesh is conforming or not 309 310 Level: advanced 311 312 Note: 313 The discrete gradient matrix G is used to analyze the subdomain edges, and it should not contain any zero entry. 314 For variable order spaces, the order should be set to zero. 315 If global is true, the rows of G should be given in global ordering for the whole dofs; 316 if false, the ordering should be global for the Nedelec field. 317 In the latter case, it should hold gid[i] < gid[j] iff geid[i] < geid[j], with gid the global orderding for all the dofs 318 and geid the one for the Nedelec field. 319 320 .seealso: `PCBDDC`, `PCBDDCSetDofsSplitting()`, `PCBDDCSetDofsSplittingLocal()`, `MATAIJ`, `PCBDDCSetDivergenceMat()` 321 @*/ 322 PetscErrorCode PCBDDCSetDiscreteGradient(PC pc, Mat G, PetscInt order, PetscInt field, PetscBool global, PetscBool conforming) 323 { 324 PetscFunctionBegin; 325 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 326 PetscValidHeaderSpecific(G, MAT_CLASSID, 2); 327 PetscValidLogicalCollectiveInt(pc, order, 3); 328 PetscValidLogicalCollectiveInt(pc, field, 4); 329 PetscValidLogicalCollectiveBool(pc, global, 5); 330 PetscValidLogicalCollectiveBool(pc, conforming, 6); 331 PetscCheckSameComm(pc, 1, G, 2); 332 PetscTryMethod(pc, "PCBDDCSetDiscreteGradient_C", (PC, Mat, PetscInt, PetscInt, PetscBool, PetscBool), (pc, G, order, field, global, conforming)); 333 PetscFunctionReturn(0); 334 } 335 336 static PetscErrorCode PCBDDCSetDivergenceMat_BDDC(PC pc, Mat divudotp, PetscBool trans, IS vl2l) 337 { 338 PC_BDDC *pcbddc = (PC_BDDC *)pc->data; 339 340 PetscFunctionBegin; 341 PetscCall(PetscObjectReference((PetscObject)divudotp)); 342 PetscCall(MatDestroy(&pcbddc->divudotp)); 343 pcbddc->divudotp = divudotp; 344 pcbddc->divudotp_trans = trans; 345 pcbddc->compute_nonetflux = PETSC_TRUE; 346 if (vl2l) { 347 PetscCall(PetscObjectReference((PetscObject)vl2l)); 348 PetscCall(ISDestroy(&pcbddc->divudotp_vl2l)); 349 pcbddc->divudotp_vl2l = vl2l; 350 } 351 PetscFunctionReturn(0); 352 } 353 354 /*@ 355 PCBDDCSetDivergenceMat - Sets the linear operator representing \int_\Omega \div {\bf u} \cdot p dx 356 357 Collective on pc 358 359 Input Parameters: 360 + pc - the preconditioning context 361 . divudotp - the matrix (must be of type `MATIS`) 362 . trans - if trans if false (resp. true), then pressures are in the test (trial) space and velocities are in the trial (test) space. 363 - vl2l - optional index set describing the local (wrt the local matrix in divudotp) to local (wrt the local matrix 364 in the preconditioning matrix) map for the velocities 365 366 Level: advanced 367 368 Notes: 369 This auxiliary matrix is used to compute quadrature weights representing the net-flux across subdomain boundaries 370 371 If vl2l is NULL, the local ordering for velocities in divudotp should match that of the preconditioning matrix 372 373 .seealso: `PCBDDC`, `PCBDDCSetDiscreteGradient()` 374 @*/ 375 PetscErrorCode PCBDDCSetDivergenceMat(PC pc, Mat divudotp, PetscBool trans, IS vl2l) 376 { 377 PetscBool ismatis; 378 379 PetscFunctionBegin; 380 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 381 PetscValidHeaderSpecific(divudotp, MAT_CLASSID, 2); 382 PetscCheckSameComm(pc, 1, divudotp, 2); 383 PetscValidLogicalCollectiveBool(pc, trans, 3); 384 if (vl2l) PetscValidHeaderSpecific(vl2l, IS_CLASSID, 4); 385 PetscCall(PetscObjectTypeCompare((PetscObject)divudotp, MATIS, &ismatis)); 386 PetscCheck(ismatis, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Divergence matrix needs to be of type MATIS"); 387 PetscTryMethod(pc, "PCBDDCSetDivergenceMat_C", (PC, Mat, PetscBool, IS), (pc, divudotp, trans, vl2l)); 388 PetscFunctionReturn(0); 389 } 390 391 static PetscErrorCode PCBDDCSetChangeOfBasisMat_BDDC(PC pc, Mat change, PetscBool interior) 392 { 393 PC_BDDC *pcbddc = (PC_BDDC *)pc->data; 394 395 PetscFunctionBegin; 396 PetscCall(PetscObjectReference((PetscObject)change)); 397 PetscCall(MatDestroy(&pcbddc->user_ChangeOfBasisMatrix)); 398 pcbddc->user_ChangeOfBasisMatrix = change; 399 pcbddc->change_interior = interior; 400 PetscFunctionReturn(0); 401 } 402 403 /*@ 404 PCBDDCSetChangeOfBasisMat - Set user defined change of basis for dofs 405 406 Collective on pc 407 408 Input Parameters: 409 + pc - the preconditioning context 410 . change - the change of basis matrix 411 - interior - whether or not the change of basis modifies interior dofs 412 413 Level: intermediate 414 415 .seealso: `PCBDDC` 416 @*/ 417 PetscErrorCode PCBDDCSetChangeOfBasisMat(PC pc, Mat change, PetscBool interior) 418 { 419 PetscFunctionBegin; 420 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 421 PetscValidHeaderSpecific(change, MAT_CLASSID, 2); 422 PetscCheckSameComm(pc, 1, change, 2); 423 if (pc->mat) { 424 PetscInt rows_c, cols_c, rows, cols; 425 PetscCall(MatGetSize(pc->mat, &rows, &cols)); 426 PetscCall(MatGetSize(change, &rows_c, &cols_c)); 427 PetscCheck(rows_c == rows, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Invalid number of rows for change of basis matrix! %" PetscInt_FMT " != %" PetscInt_FMT, rows_c, rows); 428 PetscCheck(cols_c == cols, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Invalid number of columns for change of basis matrix! %" PetscInt_FMT " != %" PetscInt_FMT, cols_c, cols); 429 PetscCall(MatGetLocalSize(pc->mat, &rows, &cols)); 430 PetscCall(MatGetLocalSize(change, &rows_c, &cols_c)); 431 PetscCheck(rows_c == rows, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Invalid number of local rows for change of basis matrix! %" PetscInt_FMT " != %" PetscInt_FMT, rows_c, rows); 432 PetscCheck(cols_c == cols, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Invalid number of local columns for change of basis matrix! %" PetscInt_FMT " != %" PetscInt_FMT, cols_c, cols); 433 } 434 PetscTryMethod(pc, "PCBDDCSetChangeOfBasisMat_C", (PC, Mat, PetscBool), (pc, change, interior)); 435 PetscFunctionReturn(0); 436 } 437 438 static PetscErrorCode PCBDDCSetPrimalVerticesIS_BDDC(PC pc, IS PrimalVertices) 439 { 440 PC_BDDC *pcbddc = (PC_BDDC *)pc->data; 441 PetscBool isequal = PETSC_FALSE; 442 443 PetscFunctionBegin; 444 PetscCall(PetscObjectReference((PetscObject)PrimalVertices)); 445 if (pcbddc->user_primal_vertices) PetscCall(ISEqual(PrimalVertices, pcbddc->user_primal_vertices, &isequal)); 446 PetscCall(ISDestroy(&pcbddc->user_primal_vertices)); 447 PetscCall(ISDestroy(&pcbddc->user_primal_vertices_local)); 448 pcbddc->user_primal_vertices = PrimalVertices; 449 if (!isequal) pcbddc->recompute_topography = PETSC_TRUE; 450 PetscFunctionReturn(0); 451 } 452 453 /*@ 454 PCBDDCSetPrimalVerticesIS - Set additional user defined primal vertices in `PCBDDC` 455 456 Collective 457 458 Input Parameters: 459 + pc - the preconditioning context 460 - PrimalVertices - index set of primal vertices in global numbering (can be empty) 461 462 Level: intermediate 463 464 Note: 465 Any process can list any global node 466 467 .seealso: `PCBDDC`, `PCBDDCGetPrimalVerticesIS()`, `PCBDDCSetPrimalVerticesLocalIS()`, `PCBDDCGetPrimalVerticesLocalIS()` 468 @*/ 469 PetscErrorCode PCBDDCSetPrimalVerticesIS(PC pc, IS PrimalVertices) 470 { 471 PetscFunctionBegin; 472 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 473 PetscValidHeaderSpecific(PrimalVertices, IS_CLASSID, 2); 474 PetscCheckSameComm(pc, 1, PrimalVertices, 2); 475 PetscTryMethod(pc, "PCBDDCSetPrimalVerticesIS_C", (PC, IS), (pc, PrimalVertices)); 476 PetscFunctionReturn(0); 477 } 478 479 static PetscErrorCode PCBDDCGetPrimalVerticesIS_BDDC(PC pc, IS *is) 480 { 481 PC_BDDC *pcbddc = (PC_BDDC *)pc->data; 482 483 PetscFunctionBegin; 484 *is = pcbddc->user_primal_vertices; 485 PetscFunctionReturn(0); 486 } 487 488 /*@ 489 PCBDDCGetPrimalVerticesIS - Get user defined primal vertices set with `PCBDDCSetPrimalVerticesIS()` 490 491 Collective 492 493 Input Parameter: 494 . pc - the preconditioning context 495 496 Output Parameter: 497 . is - index set of primal vertices in global numbering (NULL if not set) 498 499 Level: intermediate 500 501 .seealso: `PCBDDC`, `PCBDDCSetPrimalVerticesIS()`, `PCBDDCSetPrimalVerticesLocalIS()`, `PCBDDCGetPrimalVerticesLocalIS()` 502 @*/ 503 PetscErrorCode PCBDDCGetPrimalVerticesIS(PC pc, IS *is) 504 { 505 PetscFunctionBegin; 506 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 507 PetscValidPointer(is, 2); 508 PetscUseMethod(pc, "PCBDDCGetPrimalVerticesIS_C", (PC, IS *), (pc, is)); 509 PetscFunctionReturn(0); 510 } 511 512 static PetscErrorCode PCBDDCSetPrimalVerticesLocalIS_BDDC(PC pc, IS PrimalVertices) 513 { 514 PC_BDDC *pcbddc = (PC_BDDC *)pc->data; 515 PetscBool isequal = PETSC_FALSE; 516 517 PetscFunctionBegin; 518 PetscCall(PetscObjectReference((PetscObject)PrimalVertices)); 519 if (pcbddc->user_primal_vertices_local) PetscCall(ISEqual(PrimalVertices, pcbddc->user_primal_vertices_local, &isequal)); 520 PetscCall(ISDestroy(&pcbddc->user_primal_vertices)); 521 PetscCall(ISDestroy(&pcbddc->user_primal_vertices_local)); 522 pcbddc->user_primal_vertices_local = PrimalVertices; 523 if (!isequal) pcbddc->recompute_topography = PETSC_TRUE; 524 PetscFunctionReturn(0); 525 } 526 527 /*@ 528 PCBDDCSetPrimalVerticesLocalIS - Set additional user defined primal vertices in `PCBDDC` 529 530 Collective 531 532 Input Parameters: 533 + pc - the preconditioning context 534 - PrimalVertices - index set of primal vertices in local numbering (can be empty) 535 536 Level: intermediate 537 538 .seealso: `PCBDDC`, `PCBDDCSetPrimalVerticesIS()`, `PCBDDCGetPrimalVerticesIS()`, `PCBDDCGetPrimalVerticesLocalIS()` 539 @*/ 540 PetscErrorCode PCBDDCSetPrimalVerticesLocalIS(PC pc, IS PrimalVertices) 541 { 542 PetscFunctionBegin; 543 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 544 PetscValidHeaderSpecific(PrimalVertices, IS_CLASSID, 2); 545 PetscCheckSameComm(pc, 1, PrimalVertices, 2); 546 PetscTryMethod(pc, "PCBDDCSetPrimalVerticesLocalIS_C", (PC, IS), (pc, PrimalVertices)); 547 PetscFunctionReturn(0); 548 } 549 550 static PetscErrorCode PCBDDCGetPrimalVerticesLocalIS_BDDC(PC pc, IS *is) 551 { 552 PC_BDDC *pcbddc = (PC_BDDC *)pc->data; 553 554 PetscFunctionBegin; 555 *is = pcbddc->user_primal_vertices_local; 556 PetscFunctionReturn(0); 557 } 558 559 /*@ 560 PCBDDCGetPrimalVerticesLocalIS - Get user defined primal vertices set with `PCBDDCSetPrimalVerticesLocalIS()` 561 562 Collective 563 564 Input Parameter: 565 . pc - the preconditioning context 566 567 Output Parameter: 568 . is - index set of primal vertices in local numbering (NULL if not set) 569 570 Level: intermediate 571 572 .seealso: `PCBDDC`, `PCBDDCSetPrimalVerticesIS()`, `PCBDDCGetPrimalVerticesIS()`, `PCBDDCSetPrimalVerticesLocalIS()` 573 @*/ 574 PetscErrorCode PCBDDCGetPrimalVerticesLocalIS(PC pc, IS *is) 575 { 576 PetscFunctionBegin; 577 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 578 PetscValidPointer(is, 2); 579 PetscUseMethod(pc, "PCBDDCGetPrimalVerticesLocalIS_C", (PC, IS *), (pc, is)); 580 PetscFunctionReturn(0); 581 } 582 583 static PetscErrorCode PCBDDCSetCoarseningRatio_BDDC(PC pc, PetscInt k) 584 { 585 PC_BDDC *pcbddc = (PC_BDDC *)pc->data; 586 587 PetscFunctionBegin; 588 pcbddc->coarsening_ratio = k; 589 PetscFunctionReturn(0); 590 } 591 592 /*@ 593 PCBDDCSetCoarseningRatio - Set coarsening ratio used in multilevel version 594 595 Logically collective on pc 596 597 Input Parameters: 598 + pc - the preconditioning context 599 - k - coarsening ratio (H/h at the coarser level) 600 601 Options Database Key: 602 . -pc_bddc_coarsening_ratio <int> - Set coarsening ratio used in multilevel coarsening 603 604 Level: intermediate 605 606 Note: 607 Approximately k subdomains at the finer level will be aggregated into a single subdomain at the coarser level 608 609 .seealso: `PCBDDC`, `PCBDDCSetLevels()` 610 @*/ 611 PetscErrorCode PCBDDCSetCoarseningRatio(PC pc, PetscInt k) 612 { 613 PetscFunctionBegin; 614 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 615 PetscValidLogicalCollectiveInt(pc, k, 2); 616 PetscTryMethod(pc, "PCBDDCSetCoarseningRatio_C", (PC, PetscInt), (pc, k)); 617 PetscFunctionReturn(0); 618 } 619 620 /* The following functions (PCBDDCSetUseExactDirichlet PCBDDCSetLevel) are not public */ 621 static PetscErrorCode PCBDDCSetUseExactDirichlet_BDDC(PC pc, PetscBool flg) 622 { 623 PC_BDDC *pcbddc = (PC_BDDC *)pc->data; 624 625 PetscFunctionBegin; 626 pcbddc->use_exact_dirichlet_trick = flg; 627 PetscFunctionReturn(0); 628 } 629 630 PetscErrorCode PCBDDCSetUseExactDirichlet(PC pc, PetscBool flg) 631 { 632 PetscFunctionBegin; 633 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 634 PetscValidLogicalCollectiveBool(pc, flg, 2); 635 PetscTryMethod(pc, "PCBDDCSetUseExactDirichlet_C", (PC, PetscBool), (pc, flg)); 636 PetscFunctionReturn(0); 637 } 638 639 static PetscErrorCode PCBDDCSetLevel_BDDC(PC pc, PetscInt level) 640 { 641 PC_BDDC *pcbddc = (PC_BDDC *)pc->data; 642 643 PetscFunctionBegin; 644 pcbddc->current_level = level; 645 PetscFunctionReturn(0); 646 } 647 648 PetscErrorCode PCBDDCSetLevel(PC pc, PetscInt level) 649 { 650 PetscFunctionBegin; 651 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 652 PetscValidLogicalCollectiveInt(pc, level, 2); 653 PetscTryMethod(pc, "PCBDDCSetLevel_C", (PC, PetscInt), (pc, level)); 654 PetscFunctionReturn(0); 655 } 656 657 static PetscErrorCode PCBDDCSetLevels_BDDC(PC pc, PetscInt levels) 658 { 659 PC_BDDC *pcbddc = (PC_BDDC *)pc->data; 660 661 PetscFunctionBegin; 662 PetscCheck(levels < PETSC_PCBDDC_MAXLEVELS, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Maximum number of additional levels for BDDC is %d", PETSC_PCBDDC_MAXLEVELS - 1); 663 pcbddc->max_levels = levels; 664 PetscFunctionReturn(0); 665 } 666 667 /*@ 668 PCBDDCSetLevels - Sets the maximum number of additional levels allowed for multilevel `PCBDDC` 669 670 Logically collective on pc 671 672 Input Parameters: 673 + pc - the preconditioning context 674 - levels - the maximum number of levels 675 676 Options Database Key: 677 . -pc_bddc_levels <int> - Set maximum number of levels for multilevel 678 679 Level: intermediate 680 681 Note: 682 The default value is 0, that gives the classical two-levels BDDC 683 684 .seealso: `PCBDDC`, `PCBDDCSetCoarseningRatio()` 685 @*/ 686 PetscErrorCode PCBDDCSetLevels(PC pc, PetscInt levels) 687 { 688 PetscFunctionBegin; 689 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 690 PetscValidLogicalCollectiveInt(pc, levels, 2); 691 PetscTryMethod(pc, "PCBDDCSetLevels_C", (PC, PetscInt), (pc, levels)); 692 PetscFunctionReturn(0); 693 } 694 695 static PetscErrorCode PCBDDCSetDirichletBoundaries_BDDC(PC pc, IS DirichletBoundaries) 696 { 697 PC_BDDC *pcbddc = (PC_BDDC *)pc->data; 698 PetscBool isequal = PETSC_FALSE; 699 700 PetscFunctionBegin; 701 PetscCall(PetscObjectReference((PetscObject)DirichletBoundaries)); 702 if (pcbddc->DirichletBoundaries) PetscCall(ISEqual(DirichletBoundaries, pcbddc->DirichletBoundaries, &isequal)); 703 /* last user setting takes precedence -> destroy any other customization */ 704 PetscCall(ISDestroy(&pcbddc->DirichletBoundariesLocal)); 705 PetscCall(ISDestroy(&pcbddc->DirichletBoundaries)); 706 pcbddc->DirichletBoundaries = DirichletBoundaries; 707 if (!isequal) pcbddc->recompute_topography = PETSC_TRUE; 708 PetscFunctionReturn(0); 709 } 710 711 /*@ 712 PCBDDCSetDirichletBoundaries - Set IS defining Dirichlet boundaries for the global problem. 713 714 Collective 715 716 Input Parameters: 717 + pc - the preconditioning context 718 - DirichletBoundaries - parallel IS defining the Dirichlet boundaries 719 720 Level: intermediate 721 722 Note: 723 Provide the information if you used `MatZeroRows()` or `MatZeroRowsColumns()`. Any process can list any global node 724 725 .seealso: `PCBDDC`, `PCBDDCSetDirichletBoundariesLocal()`, `MatZeroRows()`, `MatZeroRowsColumns()` 726 @*/ 727 PetscErrorCode PCBDDCSetDirichletBoundaries(PC pc, IS DirichletBoundaries) 728 { 729 PetscFunctionBegin; 730 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 731 PetscValidHeaderSpecific(DirichletBoundaries, IS_CLASSID, 2); 732 PetscCheckSameComm(pc, 1, DirichletBoundaries, 2); 733 PetscTryMethod(pc, "PCBDDCSetDirichletBoundaries_C", (PC, IS), (pc, DirichletBoundaries)); 734 PetscFunctionReturn(0); 735 } 736 737 static PetscErrorCode PCBDDCSetDirichletBoundariesLocal_BDDC(PC pc, IS DirichletBoundaries) 738 { 739 PC_BDDC *pcbddc = (PC_BDDC *)pc->data; 740 PetscBool isequal = PETSC_FALSE; 741 742 PetscFunctionBegin; 743 PetscCall(PetscObjectReference((PetscObject)DirichletBoundaries)); 744 if (pcbddc->DirichletBoundariesLocal) PetscCall(ISEqual(DirichletBoundaries, pcbddc->DirichletBoundariesLocal, &isequal)); 745 /* last user setting takes precedence -> destroy any other customization */ 746 PetscCall(ISDestroy(&pcbddc->DirichletBoundariesLocal)); 747 PetscCall(ISDestroy(&pcbddc->DirichletBoundaries)); 748 pcbddc->DirichletBoundariesLocal = DirichletBoundaries; 749 if (!isequal) pcbddc->recompute_topography = PETSC_TRUE; 750 PetscFunctionReturn(0); 751 } 752 753 /*@ 754 PCBDDCSetDirichletBoundariesLocal - Set IS defining Dirichlet boundaries for the global problem in local ordering. 755 756 Collective 757 758 Input Parameters: 759 + pc - the preconditioning context 760 - DirichletBoundaries - parallel IS defining the Dirichlet boundaries (in local ordering) 761 762 Level: intermediate 763 764 .seealso: `PCBDDC`, `PCBDDCSetDirichletBoundaries()`, `MatZeroRows()`, `MatZeroRowsColumns()` 765 @*/ 766 PetscErrorCode PCBDDCSetDirichletBoundariesLocal(PC pc, IS DirichletBoundaries) 767 { 768 PetscFunctionBegin; 769 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 770 PetscValidHeaderSpecific(DirichletBoundaries, IS_CLASSID, 2); 771 PetscCheckSameComm(pc, 1, DirichletBoundaries, 2); 772 PetscTryMethod(pc, "PCBDDCSetDirichletBoundariesLocal_C", (PC, IS), (pc, DirichletBoundaries)); 773 PetscFunctionReturn(0); 774 } 775 776 static PetscErrorCode PCBDDCSetNeumannBoundaries_BDDC(PC pc, IS NeumannBoundaries) 777 { 778 PC_BDDC *pcbddc = (PC_BDDC *)pc->data; 779 PetscBool isequal = PETSC_FALSE; 780 781 PetscFunctionBegin; 782 PetscCall(PetscObjectReference((PetscObject)NeumannBoundaries)); 783 if (pcbddc->NeumannBoundaries) PetscCall(ISEqual(NeumannBoundaries, pcbddc->NeumannBoundaries, &isequal)); 784 /* last user setting takes precedence -> destroy any other customization */ 785 PetscCall(ISDestroy(&pcbddc->NeumannBoundariesLocal)); 786 PetscCall(ISDestroy(&pcbddc->NeumannBoundaries)); 787 pcbddc->NeumannBoundaries = NeumannBoundaries; 788 if (!isequal) pcbddc->recompute_topography = PETSC_TRUE; 789 PetscFunctionReturn(0); 790 } 791 792 /*@ 793 PCBDDCSetNeumannBoundaries - Set `IS` defining Neumann boundaries for the global problem. 794 795 Collective on pc 796 797 Input Parameters: 798 + pc - the preconditioning context 799 - NeumannBoundaries - parallel IS defining the Neumann boundaries 800 801 Level: intermediate 802 803 Note: 804 Any process can list any global node 805 806 .seealso: `PCBDDC`, `PCBDDCSetNeumannBoundariesLocal()` 807 @*/ 808 PetscErrorCode PCBDDCSetNeumannBoundaries(PC pc, IS NeumannBoundaries) 809 { 810 PetscFunctionBegin; 811 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 812 PetscValidHeaderSpecific(NeumannBoundaries, IS_CLASSID, 2); 813 PetscCheckSameComm(pc, 1, NeumannBoundaries, 2); 814 PetscTryMethod(pc, "PCBDDCSetNeumannBoundaries_C", (PC, IS), (pc, NeumannBoundaries)); 815 PetscFunctionReturn(0); 816 } 817 818 static PetscErrorCode PCBDDCSetNeumannBoundariesLocal_BDDC(PC pc, IS NeumannBoundaries) 819 { 820 PC_BDDC *pcbddc = (PC_BDDC *)pc->data; 821 PetscBool isequal = PETSC_FALSE; 822 823 PetscFunctionBegin; 824 PetscCall(PetscObjectReference((PetscObject)NeumannBoundaries)); 825 if (pcbddc->NeumannBoundariesLocal) PetscCall(ISEqual(NeumannBoundaries, pcbddc->NeumannBoundariesLocal, &isequal)); 826 /* last user setting takes precedence -> destroy any other customization */ 827 PetscCall(ISDestroy(&pcbddc->NeumannBoundariesLocal)); 828 PetscCall(ISDestroy(&pcbddc->NeumannBoundaries)); 829 pcbddc->NeumannBoundariesLocal = NeumannBoundaries; 830 if (!isequal) pcbddc->recompute_topography = PETSC_TRUE; 831 PetscFunctionReturn(0); 832 } 833 834 /*@ 835 PCBDDCSetNeumannBoundariesLocal - Set IS defining Neumann boundaries for the global problem in local ordering. 836 837 Collective 838 839 Input Parameters: 840 + pc - the preconditioning context 841 - NeumannBoundaries - parallel IS defining the subdomain part of Neumann boundaries (in local ordering) 842 843 Level: intermediate 844 845 .seealso: `PCBDDC`, `PCBDDCSetNeumannBoundaries()`, `PCBDDCGetDirichletBoundaries()` 846 @*/ 847 PetscErrorCode PCBDDCSetNeumannBoundariesLocal(PC pc, IS NeumannBoundaries) 848 { 849 PetscFunctionBegin; 850 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 851 PetscValidHeaderSpecific(NeumannBoundaries, IS_CLASSID, 2); 852 PetscCheckSameComm(pc, 1, NeumannBoundaries, 2); 853 PetscTryMethod(pc, "PCBDDCSetNeumannBoundariesLocal_C", (PC, IS), (pc, NeumannBoundaries)); 854 PetscFunctionReturn(0); 855 } 856 857 static PetscErrorCode PCBDDCGetDirichletBoundaries_BDDC(PC pc, IS *DirichletBoundaries) 858 { 859 PC_BDDC *pcbddc = (PC_BDDC *)pc->data; 860 861 PetscFunctionBegin; 862 *DirichletBoundaries = pcbddc->DirichletBoundaries; 863 PetscFunctionReturn(0); 864 } 865 866 /*@ 867 PCBDDCGetDirichletBoundaries - Get parallel `IS` for Dirichlet boundaries 868 869 Collective 870 871 Input Parameter: 872 . pc - the preconditioning context 873 874 Output Parameter: 875 . DirichletBoundaries - index set defining the Dirichlet boundaries 876 877 Level: intermediate 878 879 Note: 880 The `IS` returned (if any) is the same passed in earlier by the user with `PCBDDCSetDirichletBoundaries()` 881 882 .seealso: `PCBDDC`, `PCBDDCSetDirichletBoundaries()` 883 @*/ 884 PetscErrorCode PCBDDCGetDirichletBoundaries(PC pc, IS *DirichletBoundaries) 885 { 886 PetscFunctionBegin; 887 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 888 PetscUseMethod(pc, "PCBDDCGetDirichletBoundaries_C", (PC, IS *), (pc, DirichletBoundaries)); 889 PetscFunctionReturn(0); 890 } 891 892 static PetscErrorCode PCBDDCGetDirichletBoundariesLocal_BDDC(PC pc, IS *DirichletBoundaries) 893 { 894 PC_BDDC *pcbddc = (PC_BDDC *)pc->data; 895 896 PetscFunctionBegin; 897 *DirichletBoundaries = pcbddc->DirichletBoundariesLocal; 898 PetscFunctionReturn(0); 899 } 900 901 /*@ 902 PCBDDCGetDirichletBoundariesLocal - Get parallel `IS` for Dirichlet boundaries (in local ordering) 903 904 Collective 905 906 Input Parameter: 907 . pc - the preconditioning context 908 909 Output Parameter: 910 . DirichletBoundaries - index set defining the subdomain part of Dirichlet boundaries 911 912 Level: intermediate 913 914 Note: 915 The `IS` returned could be the same passed in earlier by the user (if provided with `PCBDDCSetDirichletBoundariesLocal()`) 916 or a global-to-local map of the global `IS` (if provided with `PCBDDCSetDirichletBoundaries()`). 917 In the latter case, the `IS` will be available only after `PCSetUp()`. 918 919 .seealso: `PCBDDC`, `PCBDDCGetDirichletBoundariesLocal()`, `PCBDDCGetDirichletBoundaries()`, `PCBDDCSetDirichletBoundaries()` 920 @*/ 921 PetscErrorCode PCBDDCGetDirichletBoundariesLocal(PC pc, IS *DirichletBoundaries) 922 { 923 PetscFunctionBegin; 924 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 925 PetscUseMethod(pc, "PCBDDCGetDirichletBoundariesLocal_C", (PC, IS *), (pc, DirichletBoundaries)); 926 PetscFunctionReturn(0); 927 } 928 929 static PetscErrorCode PCBDDCGetNeumannBoundaries_BDDC(PC pc, IS *NeumannBoundaries) 930 { 931 PC_BDDC *pcbddc = (PC_BDDC *)pc->data; 932 933 PetscFunctionBegin; 934 *NeumannBoundaries = pcbddc->NeumannBoundaries; 935 PetscFunctionReturn(0); 936 } 937 938 /*@ 939 PCBDDCGetNeumannBoundaries - Get parallel `IS` for Neumann boundaries 940 941 Not Collective 942 943 Input Parameter: 944 . pc - the preconditioning context 945 946 Output Parameter: 947 . NeumannBoundaries - index set defining the Neumann boundaries 948 949 Level: intermediate 950 951 Note: 952 The `IS` returned (if any) is the same passed in earlier by the user with `PCBDDCSetNeumannBoundaries()` 953 954 .seealso: `PCBDDC`, `PCBDDCSetNeumannBoundaries()`, `PCBDDCGetDirichletBoundaries()`, `PCBDDCSetDirichletBoundaries()` 955 @*/ 956 PetscErrorCode PCBDDCGetNeumannBoundaries(PC pc, IS *NeumannBoundaries) 957 { 958 PetscFunctionBegin; 959 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 960 PetscUseMethod(pc, "PCBDDCGetNeumannBoundaries_C", (PC, IS *), (pc, NeumannBoundaries)); 961 PetscFunctionReturn(0); 962 } 963 964 static PetscErrorCode PCBDDCGetNeumannBoundariesLocal_BDDC(PC pc, IS *NeumannBoundaries) 965 { 966 PC_BDDC *pcbddc = (PC_BDDC *)pc->data; 967 968 PetscFunctionBegin; 969 *NeumannBoundaries = pcbddc->NeumannBoundariesLocal; 970 PetscFunctionReturn(0); 971 } 972 973 /*@ 974 PCBDDCGetNeumannBoundariesLocal - Get parallel `IS` for Neumann boundaries (in local ordering) 975 976 Not Collective 977 978 Input Parameter: 979 . pc - the preconditioning context 980 981 Output Parameter: 982 . NeumannBoundaries - index set defining the subdomain part of Neumann boundaries 983 984 Level: intermediate 985 986 Note: 987 The `IS` returned could be the same passed in earlier by the user (if provided with `PCBDDCSetNeumannBoundariesLocal()` 988 or a global-to-local map of the global `IS` (if provided with `PCBDDCSetNeumannBoundaries()`). 989 In the latter case, the `IS` will be available after `PCSetUp()`. 990 991 .seealso: `PCBDDC``PCBDDCSetNeumannBoundaries()`, `PCBDDCSetNeumannBoundariesLocal)`, `PCBDDCGetNeumannBoundaries()` 992 @*/ 993 PetscErrorCode PCBDDCGetNeumannBoundariesLocal(PC pc, IS *NeumannBoundaries) 994 { 995 PetscFunctionBegin; 996 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 997 PetscUseMethod(pc, "PCBDDCGetNeumannBoundariesLocal_C", (PC, IS *), (pc, NeumannBoundaries)); 998 PetscFunctionReturn(0); 999 } 1000 1001 static PetscErrorCode PCBDDCSetLocalAdjacencyGraph_BDDC(PC pc, PetscInt nvtxs, const PetscInt xadj[], const PetscInt adjncy[], PetscCopyMode copymode) 1002 { 1003 PC_BDDC *pcbddc = (PC_BDDC *)pc->data; 1004 PCBDDCGraph mat_graph = pcbddc->mat_graph; 1005 PetscBool same_data = PETSC_FALSE; 1006 1007 PetscFunctionBegin; 1008 if (!nvtxs) { 1009 if (copymode == PETSC_OWN_POINTER) { 1010 PetscCall(PetscFree(xadj)); 1011 PetscCall(PetscFree(adjncy)); 1012 } 1013 PetscCall(PCBDDCGraphResetCSR(mat_graph)); 1014 PetscFunctionReturn(0); 1015 } 1016 if (mat_graph->nvtxs == nvtxs && mat_graph->freecsr) { /* we own the data */ 1017 if (mat_graph->xadj == xadj && mat_graph->adjncy == adjncy) same_data = PETSC_TRUE; 1018 if (!same_data && mat_graph->xadj[nvtxs] == xadj[nvtxs]) { 1019 PetscCall(PetscArraycmp(xadj, mat_graph->xadj, nvtxs + 1, &same_data)); 1020 if (same_data) PetscCall(PetscArraycmp(adjncy, mat_graph->adjncy, xadj[nvtxs], &same_data)); 1021 } 1022 } 1023 if (!same_data) { 1024 /* free old CSR */ 1025 PetscCall(PCBDDCGraphResetCSR(mat_graph)); 1026 /* get CSR into graph structure */ 1027 if (copymode == PETSC_COPY_VALUES) { 1028 PetscCall(PetscMalloc1(nvtxs + 1, &mat_graph->xadj)); 1029 PetscCall(PetscMalloc1(xadj[nvtxs], &mat_graph->adjncy)); 1030 PetscCall(PetscArraycpy(mat_graph->xadj, xadj, nvtxs + 1)); 1031 PetscCall(PetscArraycpy(mat_graph->adjncy, adjncy, xadj[nvtxs])); 1032 mat_graph->freecsr = PETSC_TRUE; 1033 } else if (copymode == PETSC_OWN_POINTER) { 1034 mat_graph->xadj = (PetscInt *)xadj; 1035 mat_graph->adjncy = (PetscInt *)adjncy; 1036 mat_graph->freecsr = PETSC_TRUE; 1037 } else if (copymode == PETSC_USE_POINTER) { 1038 mat_graph->xadj = (PetscInt *)xadj; 1039 mat_graph->adjncy = (PetscInt *)adjncy; 1040 mat_graph->freecsr = PETSC_FALSE; 1041 } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported copy mode %d", copymode); 1042 mat_graph->nvtxs_csr = nvtxs; 1043 pcbddc->recompute_topography = PETSC_TRUE; 1044 } 1045 PetscFunctionReturn(0); 1046 } 1047 1048 /*@ 1049 PCBDDCSetLocalAdjacencyGraph - Set adjacency structure (CSR graph) of the local degrees of freedom. 1050 1051 Not collective 1052 1053 Input Parameters: 1054 + pc - the preconditioning context. 1055 . nvtxs - number of local vertices of the graph (i.e., the number of local dofs). 1056 . xadj, adjncy - the connectivity of the dofs in CSR format. 1057 - copymode - supported modes are `PETSC_COPY_VALUES`, `PETSC_USE_POINTER` or `PETSC_OWN_POINTER`. 1058 1059 Level: intermediate 1060 1061 Note: 1062 A dof is considered connected with all local dofs if xadj[dof+1]-xadj[dof] == 1 and adjncy[xadj[dof]] is negative. 1063 1064 .seealso: `PCBDDC`, `PetscCopyMode` 1065 @*/ 1066 PetscErrorCode PCBDDCSetLocalAdjacencyGraph(PC pc, PetscInt nvtxs, const PetscInt xadj[], const PetscInt adjncy[], PetscCopyMode copymode) 1067 { 1068 void (*f)(void) = NULL; 1069 1070 PetscFunctionBegin; 1071 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1072 if (nvtxs) { 1073 PetscValidIntPointer(xadj, 3); 1074 if (xadj[nvtxs]) PetscValidIntPointer(adjncy, 4); 1075 } 1076 PetscTryMethod(pc, "PCBDDCSetLocalAdjacencyGraph_C", (PC, PetscInt, const PetscInt[], const PetscInt[], PetscCopyMode), (pc, nvtxs, xadj, adjncy, copymode)); 1077 /* free arrays if PCBDDC is not the PC type */ 1078 PetscCall(PetscObjectQueryFunction((PetscObject)pc, "PCBDDCSetLocalAdjacencyGraph_C", &f)); 1079 if (!f && copymode == PETSC_OWN_POINTER) { 1080 PetscCall(PetscFree(xadj)); 1081 PetscCall(PetscFree(adjncy)); 1082 } 1083 PetscFunctionReturn(0); 1084 } 1085 1086 static PetscErrorCode PCBDDCSetDofsSplittingLocal_BDDC(PC pc, PetscInt n_is, IS ISForDofs[]) 1087 { 1088 PC_BDDC *pcbddc = (PC_BDDC *)pc->data; 1089 PetscInt i; 1090 PetscBool isequal = PETSC_FALSE; 1091 1092 PetscFunctionBegin; 1093 if (pcbddc->n_ISForDofsLocal == n_is) { 1094 for (i = 0; i < n_is; i++) { 1095 PetscBool isequalt; 1096 PetscCall(ISEqual(ISForDofs[i], pcbddc->ISForDofsLocal[i], &isequalt)); 1097 if (!isequalt) break; 1098 } 1099 if (i == n_is) isequal = PETSC_TRUE; 1100 } 1101 for (i = 0; i < n_is; i++) PetscCall(PetscObjectReference((PetscObject)ISForDofs[i])); 1102 /* Destroy ISes if they were already set */ 1103 for (i = 0; i < pcbddc->n_ISForDofsLocal; i++) PetscCall(ISDestroy(&pcbddc->ISForDofsLocal[i])); 1104 PetscCall(PetscFree(pcbddc->ISForDofsLocal)); 1105 /* last user setting takes precedence -> destroy any other customization */ 1106 for (i = 0; i < pcbddc->n_ISForDofs; i++) PetscCall(ISDestroy(&pcbddc->ISForDofs[i])); 1107 PetscCall(PetscFree(pcbddc->ISForDofs)); 1108 pcbddc->n_ISForDofs = 0; 1109 /* allocate space then set */ 1110 if (n_is) PetscCall(PetscMalloc1(n_is, &pcbddc->ISForDofsLocal)); 1111 for (i = 0; i < n_is; i++) pcbddc->ISForDofsLocal[i] = ISForDofs[i]; 1112 pcbddc->n_ISForDofsLocal = n_is; 1113 if (n_is) pcbddc->user_provided_isfordofs = PETSC_TRUE; 1114 if (!isequal) pcbddc->recompute_topography = PETSC_TRUE; 1115 PetscFunctionReturn(0); 1116 } 1117 1118 /*@ 1119 PCBDDCSetDofsSplittingLocal - Set index sets defining fields of the local subdomain matrix 1120 1121 Collective 1122 1123 Input Parameters: 1124 + pc - the preconditioning context 1125 . n_is - number of index sets defining the fields, must be the same on all MPI ranks 1126 - ISForDofs - array of `IS` describing the fields in local ordering 1127 1128 Level: intermediate 1129 1130 Note: 1131 Not all nodes need to be listed: unlisted nodes will belong to the complement field. 1132 1133 .seealso: `PCBDDC`, `PCBDDCSetDofsSplitting()` 1134 @*/ 1135 PetscErrorCode PCBDDCSetDofsSplittingLocal(PC pc, PetscInt n_is, IS ISForDofs[]) 1136 { 1137 PetscInt i; 1138 1139 PetscFunctionBegin; 1140 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1141 PetscValidLogicalCollectiveInt(pc, n_is, 2); 1142 for (i = 0; i < n_is; i++) { 1143 PetscCheckSameComm(pc, 1, ISForDofs[i], 3); 1144 PetscValidHeaderSpecific(ISForDofs[i], IS_CLASSID, 3); 1145 } 1146 PetscTryMethod(pc, "PCBDDCSetDofsSplittingLocal_C", (PC, PetscInt, IS[]), (pc, n_is, ISForDofs)); 1147 PetscFunctionReturn(0); 1148 } 1149 1150 static PetscErrorCode PCBDDCSetDofsSplitting_BDDC(PC pc, PetscInt n_is, IS ISForDofs[]) 1151 { 1152 PC_BDDC *pcbddc = (PC_BDDC *)pc->data; 1153 PetscInt i; 1154 PetscBool isequal = PETSC_FALSE; 1155 1156 PetscFunctionBegin; 1157 if (pcbddc->n_ISForDofs == n_is) { 1158 for (i = 0; i < n_is; i++) { 1159 PetscBool isequalt; 1160 PetscCall(ISEqual(ISForDofs[i], pcbddc->ISForDofs[i], &isequalt)); 1161 if (!isequalt) break; 1162 } 1163 if (i == n_is) isequal = PETSC_TRUE; 1164 } 1165 for (i = 0; i < n_is; i++) PetscCall(PetscObjectReference((PetscObject)ISForDofs[i])); 1166 /* Destroy ISes if they were already set */ 1167 for (i = 0; i < pcbddc->n_ISForDofs; i++) PetscCall(ISDestroy(&pcbddc->ISForDofs[i])); 1168 PetscCall(PetscFree(pcbddc->ISForDofs)); 1169 /* last user setting takes precedence -> destroy any other customization */ 1170 for (i = 0; i < pcbddc->n_ISForDofsLocal; i++) PetscCall(ISDestroy(&pcbddc->ISForDofsLocal[i])); 1171 PetscCall(PetscFree(pcbddc->ISForDofsLocal)); 1172 pcbddc->n_ISForDofsLocal = 0; 1173 /* allocate space then set */ 1174 if (n_is) PetscCall(PetscMalloc1(n_is, &pcbddc->ISForDofs)); 1175 for (i = 0; i < n_is; i++) pcbddc->ISForDofs[i] = ISForDofs[i]; 1176 pcbddc->n_ISForDofs = n_is; 1177 if (n_is) pcbddc->user_provided_isfordofs = PETSC_TRUE; 1178 if (!isequal) pcbddc->recompute_topography = PETSC_TRUE; 1179 PetscFunctionReturn(0); 1180 } 1181 1182 /*@ 1183 PCBDDCSetDofsSplitting - Set index sets defining fields of the global matrix 1184 1185 Collective 1186 1187 Input Parameters: 1188 + pc - the preconditioning context 1189 . n_is - number of index sets defining the fields 1190 - ISForDofs - array of IS describing the fields in global ordering 1191 1192 Level: intermediate 1193 1194 Note: 1195 Any process can list any global node. Not all nodes need to be listed: unlisted nodes will belong to the complement field. 1196 1197 .seealso: `PCBDDC`, `PCBDDCSetDofsSplittingLocal()` 1198 @*/ 1199 PetscErrorCode PCBDDCSetDofsSplitting(PC pc, PetscInt n_is, IS ISForDofs[]) 1200 { 1201 PetscInt i; 1202 1203 PetscFunctionBegin; 1204 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1205 PetscValidLogicalCollectiveInt(pc, n_is, 2); 1206 for (i = 0; i < n_is; i++) { 1207 PetscValidHeaderSpecific(ISForDofs[i], IS_CLASSID, 3); 1208 PetscCheckSameComm(pc, 1, ISForDofs[i], 3); 1209 } 1210 PetscTryMethod(pc, "PCBDDCSetDofsSplitting_C", (PC, PetscInt, IS[]), (pc, n_is, ISForDofs)); 1211 PetscFunctionReturn(0); 1212 } 1213 1214 /* 1215 PCPreSolve_BDDC - Changes the right hand side and (if necessary) the initial 1216 guess if a transformation of basis approach has been selected. 1217 1218 Input Parameter: 1219 + pc - the preconditioner context 1220 1221 Note: 1222 The interface routine PCPreSolve() is not usually called directly by 1223 the user, but instead is called by KSPSolve(). 1224 */ 1225 static PetscErrorCode PCPreSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x) 1226 { 1227 PC_BDDC *pcbddc = (PC_BDDC *)pc->data; 1228 PC_IS *pcis = (PC_IS *)(pc->data); 1229 Vec used_vec; 1230 PetscBool iscg, save_rhs = PETSC_TRUE, benign_correction_computed; 1231 1232 PetscFunctionBegin; 1233 /* if we are working with CG, one dirichlet solve can be avoided during Krylov iterations */ 1234 if (ksp) { 1235 PetscCall(PetscObjectTypeCompareAny((PetscObject)ksp, &iscg, KSPCG, KSPGROPPCG, KSPPIPECG, KSPPIPELCG, KSPPIPECGRR, "")); 1236 if (pcbddc->benign_apply_coarse_only || pcbddc->switch_static || !iscg || pc->mat != pc->pmat) PetscCall(PCBDDCSetUseExactDirichlet(pc, PETSC_FALSE)); 1237 } 1238 if (pcbddc->benign_apply_coarse_only || pcbddc->switch_static || pc->mat != pc->pmat) PetscCall(PCBDDCSetUseExactDirichlet(pc, PETSC_FALSE)); 1239 1240 /* Creates parallel work vectors used in presolve */ 1241 if (!pcbddc->original_rhs) PetscCall(VecDuplicate(pcis->vec1_global, &pcbddc->original_rhs)); 1242 if (!pcbddc->temp_solution) PetscCall(VecDuplicate(pcis->vec1_global, &pcbddc->temp_solution)); 1243 1244 pcbddc->temp_solution_used = PETSC_FALSE; 1245 if (x) { 1246 PetscCall(PetscObjectReference((PetscObject)x)); 1247 used_vec = x; 1248 } else { /* it can only happen when calling PCBDDCMatFETIDPGetRHS */ 1249 PetscCall(PetscObjectReference((PetscObject)pcbddc->temp_solution)); 1250 used_vec = pcbddc->temp_solution; 1251 PetscCall(VecSet(used_vec, 0.0)); 1252 pcbddc->temp_solution_used = PETSC_TRUE; 1253 PetscCall(VecCopy(rhs, pcbddc->original_rhs)); 1254 save_rhs = PETSC_FALSE; 1255 pcbddc->eliminate_dirdofs = PETSC_TRUE; 1256 } 1257 1258 /* hack into ksp data structure since PCPreSolve comes earlier than setting to zero the guess in src/ksp/ksp/interface/itfunc.c */ 1259 if (ksp) { 1260 /* store the flag for the initial guess since it will be restored back during PCPostSolve_BDDC */ 1261 PetscCall(KSPGetInitialGuessNonzero(ksp, &pcbddc->ksp_guess_nonzero)); 1262 if (!pcbddc->ksp_guess_nonzero) PetscCall(VecSet(used_vec, 0.0)); 1263 } 1264 1265 pcbddc->rhs_change = PETSC_FALSE; 1266 /* Take into account zeroed rows -> change rhs and store solution removed */ 1267 if (rhs && pcbddc->eliminate_dirdofs) { 1268 IS dirIS = NULL; 1269 1270 /* DirichletBoundariesLocal may not be consistent among neighbours; gets a dirichlet dofs IS from graph (may be cached) */ 1271 PetscCall(PCBDDCGraphGetDirichletDofs(pcbddc->mat_graph, &dirIS)); 1272 if (dirIS) { 1273 Mat_IS *matis = (Mat_IS *)pc->pmat->data; 1274 PetscInt dirsize, i, *is_indices; 1275 PetscScalar *array_x; 1276 const PetscScalar *array_diagonal; 1277 1278 PetscCall(MatGetDiagonal(pc->pmat, pcis->vec1_global)); 1279 PetscCall(VecPointwiseDivide(pcis->vec1_global, rhs, pcis->vec1_global)); 1280 PetscCall(VecScatterBegin(matis->rctx, pcis->vec1_global, pcis->vec2_N, INSERT_VALUES, SCATTER_FORWARD)); 1281 PetscCall(VecScatterEnd(matis->rctx, pcis->vec1_global, pcis->vec2_N, INSERT_VALUES, SCATTER_FORWARD)); 1282 PetscCall(VecScatterBegin(matis->rctx, used_vec, pcis->vec1_N, INSERT_VALUES, SCATTER_FORWARD)); 1283 PetscCall(VecScatterEnd(matis->rctx, used_vec, pcis->vec1_N, INSERT_VALUES, SCATTER_FORWARD)); 1284 PetscCall(ISGetLocalSize(dirIS, &dirsize)); 1285 PetscCall(VecGetArray(pcis->vec1_N, &array_x)); 1286 PetscCall(VecGetArrayRead(pcis->vec2_N, &array_diagonal)); 1287 PetscCall(ISGetIndices(dirIS, (const PetscInt **)&is_indices)); 1288 for (i = 0; i < dirsize; i++) array_x[is_indices[i]] = array_diagonal[is_indices[i]]; 1289 PetscCall(ISRestoreIndices(dirIS, (const PetscInt **)&is_indices)); 1290 PetscCall(VecRestoreArrayRead(pcis->vec2_N, &array_diagonal)); 1291 PetscCall(VecRestoreArray(pcis->vec1_N, &array_x)); 1292 PetscCall(VecScatterBegin(matis->rctx, pcis->vec1_N, used_vec, INSERT_VALUES, SCATTER_REVERSE)); 1293 PetscCall(VecScatterEnd(matis->rctx, pcis->vec1_N, used_vec, INSERT_VALUES, SCATTER_REVERSE)); 1294 pcbddc->rhs_change = PETSC_TRUE; 1295 PetscCall(ISDestroy(&dirIS)); 1296 } 1297 } 1298 1299 /* remove the computed solution or the initial guess from the rhs */ 1300 if (pcbddc->rhs_change || (ksp && pcbddc->ksp_guess_nonzero)) { 1301 /* save the original rhs */ 1302 if (save_rhs) { 1303 PetscCall(VecSwap(rhs, pcbddc->original_rhs)); 1304 save_rhs = PETSC_FALSE; 1305 } 1306 pcbddc->rhs_change = PETSC_TRUE; 1307 PetscCall(VecScale(used_vec, -1.0)); 1308 PetscCall(MatMultAdd(pc->mat, used_vec, pcbddc->original_rhs, rhs)); 1309 PetscCall(VecScale(used_vec, -1.0)); 1310 PetscCall(VecCopy(used_vec, pcbddc->temp_solution)); 1311 pcbddc->temp_solution_used = PETSC_TRUE; 1312 if (ksp) PetscCall(KSPSetInitialGuessNonzero(ksp, PETSC_FALSE)); 1313 } 1314 PetscCall(VecDestroy(&used_vec)); 1315 1316 /* compute initial vector in benign space if needed 1317 and remove non-benign solution from the rhs */ 1318 benign_correction_computed = PETSC_FALSE; 1319 if (rhs && pcbddc->benign_compute_correction && (pcbddc->benign_have_null || pcbddc->benign_apply_coarse_only)) { 1320 /* compute u^*_h using ideas similar to those in Xuemin Tu's PhD thesis (see Section 4.8.1) 1321 Recursively apply BDDC in the multilevel case */ 1322 if (!pcbddc->benign_vec) PetscCall(VecDuplicate(rhs, &pcbddc->benign_vec)); 1323 /* keep applying coarse solver unless we no longer have benign subdomains */ 1324 pcbddc->benign_apply_coarse_only = pcbddc->benign_have_null ? PETSC_TRUE : PETSC_FALSE; 1325 if (!pcbddc->benign_skip_correction) { 1326 PetscCall(PCApply_BDDC(pc, rhs, pcbddc->benign_vec)); 1327 benign_correction_computed = PETSC_TRUE; 1328 if (pcbddc->temp_solution_used) PetscCall(VecAXPY(pcbddc->temp_solution, 1.0, pcbddc->benign_vec)); 1329 PetscCall(VecScale(pcbddc->benign_vec, -1.0)); 1330 /* store the original rhs if not done earlier */ 1331 if (save_rhs) PetscCall(VecSwap(rhs, pcbddc->original_rhs)); 1332 if (pcbddc->rhs_change) { 1333 PetscCall(MatMultAdd(pc->mat, pcbddc->benign_vec, rhs, rhs)); 1334 } else { 1335 PetscCall(MatMultAdd(pc->mat, pcbddc->benign_vec, pcbddc->original_rhs, rhs)); 1336 } 1337 pcbddc->rhs_change = PETSC_TRUE; 1338 } 1339 pcbddc->benign_apply_coarse_only = PETSC_FALSE; 1340 } else { 1341 PetscCall(VecDestroy(&pcbddc->benign_vec)); 1342 } 1343 1344 /* dbg output */ 1345 if (pcbddc->dbg_flag && benign_correction_computed) { 1346 Vec v; 1347 1348 PetscCall(VecDuplicate(pcis->vec1_global, &v)); 1349 if (pcbddc->ChangeOfBasisMatrix) { 1350 PetscCall(MatMultTranspose(pcbddc->ChangeOfBasisMatrix, rhs, v)); 1351 } else { 1352 PetscCall(VecCopy(rhs, v)); 1353 } 1354 PetscCall(PCBDDCBenignGetOrSetP0(pc, v, PETSC_TRUE)); 1355 PetscCall(PetscViewerASCIIPrintf(pcbddc->dbg_viewer, "LEVEL %" PetscInt_FMT ": is the correction benign?\n", pcbddc->current_level)); 1356 PetscCall(PetscScalarView(pcbddc->benign_n, pcbddc->benign_p0, pcbddc->dbg_viewer)); 1357 PetscCall(PetscViewerFlush(pcbddc->dbg_viewer)); 1358 PetscCall(VecDestroy(&v)); 1359 } 1360 1361 /* set initial guess if using PCG */ 1362 pcbddc->exact_dirichlet_trick_app = PETSC_FALSE; 1363 if (x && pcbddc->use_exact_dirichlet_trick) { 1364 PetscCall(VecSet(x, 0.0)); 1365 if (pcbddc->ChangeOfBasisMatrix && pcbddc->change_interior) { 1366 if (benign_correction_computed) { /* we have already saved the changed rhs */ 1367 PetscCall(VecLockReadPop(pcis->vec1_global)); 1368 } else { 1369 PetscCall(MatMultTranspose(pcbddc->ChangeOfBasisMatrix, rhs, pcis->vec1_global)); 1370 } 1371 PetscCall(VecScatterBegin(pcis->global_to_D, pcis->vec1_global, pcis->vec1_D, INSERT_VALUES, SCATTER_FORWARD)); 1372 PetscCall(VecScatterEnd(pcis->global_to_D, pcis->vec1_global, pcis->vec1_D, INSERT_VALUES, SCATTER_FORWARD)); 1373 } else { 1374 PetscCall(VecScatterBegin(pcis->global_to_D, rhs, pcis->vec1_D, INSERT_VALUES, SCATTER_FORWARD)); 1375 PetscCall(VecScatterEnd(pcis->global_to_D, rhs, pcis->vec1_D, INSERT_VALUES, SCATTER_FORWARD)); 1376 } 1377 PetscCall(PetscLogEventBegin(PC_BDDC_Solves[pcbddc->current_level][0], pc, 0, 0, 0)); 1378 PetscCall(KSPSolve(pcbddc->ksp_D, pcis->vec1_D, pcis->vec2_D)); 1379 PetscCall(PetscLogEventEnd(PC_BDDC_Solves[pcbddc->current_level][0], pc, 0, 0, 0)); 1380 PetscCall(KSPCheckSolve(pcbddc->ksp_D, pc, pcis->vec2_D)); 1381 if (pcbddc->ChangeOfBasisMatrix && pcbddc->change_interior) { 1382 PetscCall(VecSet(pcis->vec1_global, 0.)); 1383 PetscCall(VecScatterBegin(pcis->global_to_D, pcis->vec2_D, pcis->vec1_global, INSERT_VALUES, SCATTER_REVERSE)); 1384 PetscCall(VecScatterEnd(pcis->global_to_D, pcis->vec2_D, pcis->vec1_global, INSERT_VALUES, SCATTER_REVERSE)); 1385 PetscCall(MatMult(pcbddc->ChangeOfBasisMatrix, pcis->vec1_global, x)); 1386 } else { 1387 PetscCall(VecScatterBegin(pcis->global_to_D, pcis->vec2_D, x, INSERT_VALUES, SCATTER_REVERSE)); 1388 PetscCall(VecScatterEnd(pcis->global_to_D, pcis->vec2_D, x, INSERT_VALUES, SCATTER_REVERSE)); 1389 } 1390 if (ksp) PetscCall(KSPSetInitialGuessNonzero(ksp, PETSC_TRUE)); 1391 pcbddc->exact_dirichlet_trick_app = PETSC_TRUE; 1392 } else if (pcbddc->ChangeOfBasisMatrix && pcbddc->change_interior && benign_correction_computed && pcbddc->use_exact_dirichlet_trick) { 1393 PetscCall(VecLockReadPop(pcis->vec1_global)); 1394 } 1395 PetscFunctionReturn(0); 1396 } 1397 1398 /* 1399 PCPostSolve_BDDC - Changes the computed solution if a transformation of basis 1400 approach has been selected. Also, restores rhs to its original state. 1401 1402 Input Parameter: 1403 + pc - the preconditioner context 1404 1405 Application Interface Routine: PCPostSolve() 1406 1407 Note: 1408 The interface routine PCPostSolve() is not usually called directly by 1409 the user, but instead is called by KSPSolve(). 1410 */ 1411 static PetscErrorCode PCPostSolve_BDDC(PC pc, KSP ksp, Vec rhs, Vec x) 1412 { 1413 PC_BDDC *pcbddc = (PC_BDDC *)pc->data; 1414 1415 PetscFunctionBegin; 1416 /* add solution removed in presolve */ 1417 if (x && pcbddc->rhs_change) { 1418 if (pcbddc->temp_solution_used) { 1419 PetscCall(VecAXPY(x, 1.0, pcbddc->temp_solution)); 1420 } else if (pcbddc->benign_compute_correction && pcbddc->benign_vec) { 1421 PetscCall(VecAXPY(x, -1.0, pcbddc->benign_vec)); 1422 } 1423 /* restore to original state (not for FETI-DP) */ 1424 if (ksp) pcbddc->temp_solution_used = PETSC_FALSE; 1425 } 1426 1427 /* restore rhs to its original state (not needed for FETI-DP) */ 1428 if (rhs && pcbddc->rhs_change) { 1429 PetscCall(VecSwap(rhs, pcbddc->original_rhs)); 1430 pcbddc->rhs_change = PETSC_FALSE; 1431 } 1432 /* restore ksp guess state */ 1433 if (ksp) { 1434 PetscCall(KSPSetInitialGuessNonzero(ksp, pcbddc->ksp_guess_nonzero)); 1435 /* reset flag for exact dirichlet trick */ 1436 pcbddc->exact_dirichlet_trick_app = PETSC_FALSE; 1437 } 1438 PetscFunctionReturn(0); 1439 } 1440 1441 /* 1442 PCSetUp_BDDC - Prepares for the use of the BDDC preconditioner 1443 by setting data structures and options. 1444 1445 Input Parameter: 1446 + pc - the preconditioner context 1447 1448 Application Interface Routine: PCSetUp() 1449 1450 Note: 1451 The interface routine PCSetUp() is not usually called directly by 1452 the user, but instead is called by PCApply() if necessary. 1453 */ 1454 PetscErrorCode PCSetUp_BDDC(PC pc) 1455 { 1456 PC_BDDC *pcbddc = (PC_BDDC *)pc->data; 1457 PCBDDCSubSchurs sub_schurs; 1458 Mat_IS *matis; 1459 MatNullSpace nearnullspace; 1460 Mat lA; 1461 IS lP, zerodiag = NULL; 1462 PetscInt nrows, ncols; 1463 PetscMPIInt size; 1464 PetscBool computesubschurs; 1465 PetscBool computeconstraintsmatrix; 1466 PetscBool new_nearnullspace_provided, ismatis, rl; 1467 PetscBool isset, issym, isspd; 1468 1469 PetscFunctionBegin; 1470 PetscCall(PetscObjectTypeCompare((PetscObject)pc->pmat, MATIS, &ismatis)); 1471 PetscCheck(ismatis, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "PCBDDC preconditioner requires matrix of type MATIS"); 1472 PetscCall(MatGetSize(pc->pmat, &nrows, &ncols)); 1473 PetscCheck(nrows == ncols, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "PCBDDC preconditioner requires a square preconditioning matrix"); 1474 PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size)); 1475 1476 matis = (Mat_IS *)pc->pmat->data; 1477 /* the following lines of code should be replaced by a better logic between PCIS, PCNN, PCBDDC and other future nonoverlapping preconditioners */ 1478 /* For BDDC we need to define a local "Neumann" problem different to that defined in PCISSetup 1479 Also, BDDC builds its own KSP for the Dirichlet problem */ 1480 rl = pcbddc->recompute_topography; 1481 if (!pc->setupcalled || pc->flag == DIFFERENT_NONZERO_PATTERN) rl = PETSC_TRUE; 1482 PetscCall(MPIU_Allreduce(&rl, &pcbddc->recompute_topography, 1, MPIU_BOOL, MPI_LOR, PetscObjectComm((PetscObject)pc))); 1483 if (pcbddc->recompute_topography) { 1484 pcbddc->graphanalyzed = PETSC_FALSE; 1485 computeconstraintsmatrix = PETSC_TRUE; 1486 } else { 1487 computeconstraintsmatrix = PETSC_FALSE; 1488 } 1489 1490 /* check parameters' compatibility */ 1491 if (!pcbddc->use_deluxe_scaling) pcbddc->deluxe_zerorows = PETSC_FALSE; 1492 pcbddc->adaptive_selection = (PetscBool)(pcbddc->adaptive_threshold[0] != 0.0 || pcbddc->adaptive_threshold[1] != 0.0); 1493 pcbddc->use_deluxe_scaling = (PetscBool)(pcbddc->use_deluxe_scaling && size > 1); 1494 pcbddc->adaptive_selection = (PetscBool)(pcbddc->adaptive_selection && size > 1); 1495 pcbddc->adaptive_userdefined = (PetscBool)(pcbddc->adaptive_selection && pcbddc->adaptive_userdefined); 1496 if (pcbddc->adaptive_selection) pcbddc->use_faces = PETSC_TRUE; 1497 1498 computesubschurs = (PetscBool)(pcbddc->adaptive_selection || pcbddc->use_deluxe_scaling); 1499 1500 /* activate all connected components if the netflux has been requested */ 1501 if (pcbddc->compute_nonetflux) { 1502 pcbddc->use_vertices = PETSC_TRUE; 1503 pcbddc->use_edges = PETSC_TRUE; 1504 pcbddc->use_faces = PETSC_TRUE; 1505 } 1506 1507 /* Get stdout for dbg */ 1508 if (pcbddc->dbg_flag) { 1509 if (!pcbddc->dbg_viewer) pcbddc->dbg_viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)pc)); 1510 PetscCall(PetscViewerASCIIPushSynchronized(pcbddc->dbg_viewer)); 1511 PetscCall(PetscViewerASCIIAddTab(pcbddc->dbg_viewer, 2 * pcbddc->current_level)); 1512 } 1513 1514 /* process topology information */ 1515 PetscCall(PetscLogEventBegin(PC_BDDC_Topology[pcbddc->current_level], pc, 0, 0, 0)); 1516 if (pcbddc->recompute_topography) { 1517 PetscCall(PCBDDCComputeLocalTopologyInfo(pc)); 1518 if (pcbddc->discretegradient) PetscCall(PCBDDCNedelecSupport(pc)); 1519 } 1520 if (pcbddc->corner_selected) pcbddc->use_vertices = PETSC_TRUE; 1521 1522 /* change basis if requested by the user */ 1523 if (pcbddc->user_ChangeOfBasisMatrix) { 1524 /* use_change_of_basis flag is used to automatically compute a change of basis from constraints */ 1525 pcbddc->use_change_of_basis = PETSC_FALSE; 1526 PetscCall(PCBDDCComputeLocalMatrix(pc, pcbddc->user_ChangeOfBasisMatrix)); 1527 } else { 1528 PetscCall(MatDestroy(&pcbddc->local_mat)); 1529 PetscCall(PetscObjectReference((PetscObject)matis->A)); 1530 pcbddc->local_mat = matis->A; 1531 } 1532 1533 /* 1534 Compute change of basis on local pressures (aka zerodiag dofs) with the benign trick 1535 This should come earlier then PCISSetUp for extracting the correct subdomain matrices 1536 */ 1537 PetscCall(PCBDDCBenignShellMat(pc, PETSC_TRUE)); 1538 if (pcbddc->benign_saddle_point) { 1539 PC_IS *pcis = (PC_IS *)pc->data; 1540 1541 if (pcbddc->user_ChangeOfBasisMatrix || pcbddc->use_change_of_basis || !computesubschurs) pcbddc->benign_change_explicit = PETSC_TRUE; 1542 /* detect local saddle point and change the basis in pcbddc->local_mat */ 1543 PetscCall(PCBDDCBenignDetectSaddlePoint(pc, (PetscBool)(!pcbddc->recompute_topography), &zerodiag)); 1544 /* pop B0 mat from local mat */ 1545 PetscCall(PCBDDCBenignPopOrPushB0(pc, PETSC_TRUE)); 1546 /* give pcis a hint to not reuse submatrices during PCISCreate */ 1547 if (pc->flag == SAME_NONZERO_PATTERN && pcis->reusesubmatrices == PETSC_TRUE) { 1548 if (pcbddc->benign_n && (pcbddc->benign_change_explicit || pcbddc->dbg_flag)) { 1549 pcis->reusesubmatrices = PETSC_FALSE; 1550 } else { 1551 pcis->reusesubmatrices = PETSC_TRUE; 1552 } 1553 } else { 1554 pcis->reusesubmatrices = PETSC_FALSE; 1555 } 1556 } 1557 1558 /* propagate relevant information */ 1559 PetscCall(MatIsSymmetricKnown(matis->A, &isset, &issym)); 1560 if (isset) PetscCall(MatSetOption(pcbddc->local_mat, MAT_SYMMETRIC, issym)); 1561 PetscCall(MatIsSPDKnown(matis->A, &isset, &isspd)); 1562 if (isset) PetscCall(MatSetOption(pcbddc->local_mat, MAT_SPD, isspd)); 1563 1564 /* Set up all the "iterative substructuring" common block without computing solvers */ 1565 { 1566 Mat temp_mat; 1567 1568 temp_mat = matis->A; 1569 matis->A = pcbddc->local_mat; 1570 PetscCall(PCISSetUp(pc, PETSC_TRUE, PETSC_FALSE)); 1571 pcbddc->local_mat = matis->A; 1572 matis->A = temp_mat; 1573 } 1574 1575 /* Analyze interface */ 1576 if (!pcbddc->graphanalyzed) { 1577 PetscCall(PCBDDCAnalyzeInterface(pc)); 1578 computeconstraintsmatrix = PETSC_TRUE; 1579 if (pcbddc->adaptive_selection && !pcbddc->use_deluxe_scaling && !pcbddc->mat_graph->twodim) { 1580 SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Cannot compute the adaptive primal space for a problem with 3D edges without deluxe scaling"); 1581 } 1582 if (pcbddc->compute_nonetflux) { 1583 MatNullSpace nnfnnsp; 1584 1585 PetscCheck(pcbddc->divudotp, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Missing divudotp operator"); 1586 PetscCall(PCBDDCComputeNoNetFlux(pc->pmat, pcbddc->divudotp, pcbddc->divudotp_trans, pcbddc->divudotp_vl2l, pcbddc->mat_graph, &nnfnnsp)); 1587 /* TODO what if a nearnullspace is already attached? */ 1588 if (nnfnnsp) { 1589 PetscCall(MatSetNearNullSpace(pc->pmat, nnfnnsp)); 1590 PetscCall(MatNullSpaceDestroy(&nnfnnsp)); 1591 } 1592 } 1593 } 1594 PetscCall(PetscLogEventEnd(PC_BDDC_Topology[pcbddc->current_level], pc, 0, 0, 0)); 1595 1596 /* check existence of a divergence free extension, i.e. 1597 b(v_I,p_0) = 0 for all v_I (raise error if not). 1598 Also, check that PCBDDCBenignGetOrSetP0 works */ 1599 if (pcbddc->benign_saddle_point && pcbddc->dbg_flag > 1) PetscCall(PCBDDCBenignCheck(pc, zerodiag)); 1600 PetscCall(ISDestroy(&zerodiag)); 1601 1602 /* Setup local dirichlet solver ksp_D and sub_schurs solvers */ 1603 if (computesubschurs && pcbddc->recompute_topography) PetscCall(PCBDDCInitSubSchurs(pc)); 1604 /* SetUp Scaling operator (scaling matrices could be needed in SubSchursSetUp)*/ 1605 if (!pcbddc->use_deluxe_scaling) PetscCall(PCBDDCScalingSetUp(pc)); 1606 1607 /* finish setup solvers and do adaptive selection of constraints */ 1608 sub_schurs = pcbddc->sub_schurs; 1609 if (sub_schurs && sub_schurs->schur_explicit) { 1610 if (computesubschurs) PetscCall(PCBDDCSetUpSubSchurs(pc)); 1611 PetscCall(PCBDDCSetUpLocalSolvers(pc, PETSC_TRUE, PETSC_FALSE)); 1612 } else { 1613 PetscCall(PCBDDCSetUpLocalSolvers(pc, PETSC_TRUE, PETSC_FALSE)); 1614 if (computesubschurs) PetscCall(PCBDDCSetUpSubSchurs(pc)); 1615 } 1616 if (pcbddc->adaptive_selection) { 1617 PetscCall(PCBDDCAdaptiveSelection(pc)); 1618 computeconstraintsmatrix = PETSC_TRUE; 1619 } 1620 1621 /* infer if NullSpace object attached to Mat via MatSetNearNullSpace has changed */ 1622 new_nearnullspace_provided = PETSC_FALSE; 1623 PetscCall(MatGetNearNullSpace(pc->pmat, &nearnullspace)); 1624 if (pcbddc->onearnullspace) { /* already used nearnullspace */ 1625 if (!nearnullspace) { /* near null space attached to mat has been destroyed */ 1626 new_nearnullspace_provided = PETSC_TRUE; 1627 } else { 1628 /* determine if the two nullspaces are different (should be lightweight) */ 1629 if (nearnullspace != pcbddc->onearnullspace) { 1630 new_nearnullspace_provided = PETSC_TRUE; 1631 } else { /* maybe the user has changed the content of the nearnullspace so check vectors ObjectStateId */ 1632 PetscInt i; 1633 const Vec *nearnullvecs; 1634 PetscObjectState state; 1635 PetscInt nnsp_size; 1636 PetscCall(MatNullSpaceGetVecs(nearnullspace, NULL, &nnsp_size, &nearnullvecs)); 1637 for (i = 0; i < nnsp_size; i++) { 1638 PetscCall(PetscObjectStateGet((PetscObject)nearnullvecs[i], &state)); 1639 if (pcbddc->onearnullvecs_state[i] != state) { 1640 new_nearnullspace_provided = PETSC_TRUE; 1641 break; 1642 } 1643 } 1644 } 1645 } 1646 } else { 1647 if (!nearnullspace) { /* both nearnullspaces are null */ 1648 new_nearnullspace_provided = PETSC_FALSE; 1649 } else { /* nearnullspace attached later */ 1650 new_nearnullspace_provided = PETSC_TRUE; 1651 } 1652 } 1653 1654 /* Setup constraints and related work vectors */ 1655 /* reset primal space flags */ 1656 PetscCall(PetscLogEventBegin(PC_BDDC_LocalWork[pcbddc->current_level], pc, 0, 0, 0)); 1657 pcbddc->new_primal_space = PETSC_FALSE; 1658 pcbddc->new_primal_space_local = PETSC_FALSE; 1659 if (computeconstraintsmatrix || new_nearnullspace_provided) { 1660 /* It also sets the primal space flags */ 1661 PetscCall(PCBDDCConstraintsSetUp(pc)); 1662 } 1663 /* Allocate needed local vectors (which depends on quantities defined during ConstraintsSetUp) */ 1664 PetscCall(PCBDDCSetUpLocalWorkVectors(pc)); 1665 1666 if (pcbddc->use_change_of_basis) { 1667 PC_IS *pcis = (PC_IS *)(pc->data); 1668 1669 PetscCall(PCBDDCComputeLocalMatrix(pc, pcbddc->ChangeOfBasisMatrix)); 1670 if (pcbddc->benign_change) { 1671 PetscCall(MatDestroy(&pcbddc->benign_B0)); 1672 /* pop B0 from pcbddc->local_mat */ 1673 PetscCall(PCBDDCBenignPopOrPushB0(pc, PETSC_TRUE)); 1674 } 1675 /* get submatrices */ 1676 PetscCall(MatDestroy(&pcis->A_IB)); 1677 PetscCall(MatDestroy(&pcis->A_BI)); 1678 PetscCall(MatDestroy(&pcis->A_BB)); 1679 PetscCall(MatCreateSubMatrix(pcbddc->local_mat, pcis->is_B_local, pcis->is_B_local, MAT_INITIAL_MATRIX, &pcis->A_BB)); 1680 PetscCall(MatCreateSubMatrix(pcbddc->local_mat, pcis->is_I_local, pcis->is_B_local, MAT_INITIAL_MATRIX, &pcis->A_IB)); 1681 PetscCall(MatCreateSubMatrix(pcbddc->local_mat, pcis->is_B_local, pcis->is_I_local, MAT_INITIAL_MATRIX, &pcis->A_BI)); 1682 /* set flag in pcis to not reuse submatrices during PCISCreate */ 1683 pcis->reusesubmatrices = PETSC_FALSE; 1684 } else if (!pcbddc->user_ChangeOfBasisMatrix && !pcbddc->benign_change) { 1685 PetscCall(MatDestroy(&pcbddc->local_mat)); 1686 PetscCall(PetscObjectReference((PetscObject)matis->A)); 1687 pcbddc->local_mat = matis->A; 1688 } 1689 1690 /* interface pressure block row for B_C */ 1691 PetscCall(PetscObjectQuery((PetscObject)pc, "__KSPFETIDP_lP", (PetscObject *)&lP)); 1692 PetscCall(PetscObjectQuery((PetscObject)pc, "__KSPFETIDP_lA", (PetscObject *)&lA)); 1693 if (lA && lP) { 1694 PC_IS *pcis = (PC_IS *)pc->data; 1695 Mat B_BI, B_BB, Bt_BI, Bt_BB; 1696 PetscBool issym; 1697 1698 PetscCall(MatIsSymmetric(lA, PETSC_SMALL, &issym)); 1699 if (issym) { 1700 PetscCall(MatCreateSubMatrix(lA, lP, pcis->is_I_local, MAT_INITIAL_MATRIX, &B_BI)); 1701 PetscCall(MatCreateSubMatrix(lA, lP, pcis->is_B_local, MAT_INITIAL_MATRIX, &B_BB)); 1702 PetscCall(MatCreateTranspose(B_BI, &Bt_BI)); 1703 PetscCall(MatCreateTranspose(B_BB, &Bt_BB)); 1704 } else { 1705 PetscCall(MatCreateSubMatrix(lA, lP, pcis->is_I_local, MAT_INITIAL_MATRIX, &B_BI)); 1706 PetscCall(MatCreateSubMatrix(lA, lP, pcis->is_B_local, MAT_INITIAL_MATRIX, &B_BB)); 1707 PetscCall(MatCreateSubMatrix(lA, pcis->is_I_local, lP, MAT_INITIAL_MATRIX, &Bt_BI)); 1708 PetscCall(MatCreateSubMatrix(lA, pcis->is_B_local, lP, MAT_INITIAL_MATRIX, &Bt_BB)); 1709 } 1710 PetscCall(PetscObjectCompose((PetscObject)pc, "__KSPFETIDP_B_BI", (PetscObject)B_BI)); 1711 PetscCall(PetscObjectCompose((PetscObject)pc, "__KSPFETIDP_B_BB", (PetscObject)B_BB)); 1712 PetscCall(PetscObjectCompose((PetscObject)pc, "__KSPFETIDP_Bt_BI", (PetscObject)Bt_BI)); 1713 PetscCall(PetscObjectCompose((PetscObject)pc, "__KSPFETIDP_Bt_BB", (PetscObject)Bt_BB)); 1714 PetscCall(MatDestroy(&B_BI)); 1715 PetscCall(MatDestroy(&B_BB)); 1716 PetscCall(MatDestroy(&Bt_BI)); 1717 PetscCall(MatDestroy(&Bt_BB)); 1718 } 1719 PetscCall(PetscLogEventEnd(PC_BDDC_LocalWork[pcbddc->current_level], pc, 0, 0, 0)); 1720 1721 /* SetUp coarse and local Neumann solvers */ 1722 PetscCall(PCBDDCSetUpSolvers(pc)); 1723 /* SetUp Scaling operator */ 1724 if (pcbddc->use_deluxe_scaling) PetscCall(PCBDDCScalingSetUp(pc)); 1725 1726 /* mark topography as done */ 1727 pcbddc->recompute_topography = PETSC_FALSE; 1728 1729 /* wrap pcis->A_IB and pcis->A_BI if we did not change explicitly the variables on the pressures */ 1730 PetscCall(PCBDDCBenignShellMat(pc, PETSC_FALSE)); 1731 1732 if (pcbddc->dbg_flag) { 1733 PetscCall(PetscViewerASCIISubtractTab(pcbddc->dbg_viewer, 2 * pcbddc->current_level)); 1734 PetscCall(PetscViewerASCIIPopSynchronized(pcbddc->dbg_viewer)); 1735 } 1736 PetscFunctionReturn(0); 1737 } 1738 1739 /* 1740 PCApply_BDDC - Applies the BDDC operator to a vector. 1741 1742 Input Parameters: 1743 + pc - the preconditioner context 1744 - r - input vector (global) 1745 1746 Output Parameter: 1747 . z - output vector (global) 1748 1749 Application Interface Routine: PCApply() 1750 */ 1751 PetscErrorCode PCApply_BDDC(PC pc, Vec r, Vec z) 1752 { 1753 PC_IS *pcis = (PC_IS *)(pc->data); 1754 PC_BDDC *pcbddc = (PC_BDDC *)(pc->data); 1755 Mat lA = NULL; 1756 PetscInt n_B = pcis->n_B, n_D = pcis->n - n_B; 1757 const PetscScalar one = 1.0; 1758 const PetscScalar m_one = -1.0; 1759 const PetscScalar zero = 0.0; 1760 /* This code is similar to that provided in nn.c for PCNN 1761 NN interface preconditioner changed to BDDC 1762 Added support for M_3 preconditioner in the reference article (code is active if pcbddc->switch_static == PETSC_TRUE) */ 1763 1764 PetscFunctionBegin; 1765 PetscCall(PetscCitationsRegister(citation, &cited)); 1766 if (pcbddc->switch_static) PetscCall(MatISGetLocalMat(pc->useAmat ? pc->mat : pc->pmat, &lA)); 1767 1768 if (pcbddc->ChangeOfBasisMatrix) { 1769 Vec swap; 1770 1771 PetscCall(MatMultTranspose(pcbddc->ChangeOfBasisMatrix, r, pcbddc->work_change)); 1772 swap = pcbddc->work_change; 1773 pcbddc->work_change = r; 1774 r = swap; 1775 /* save rhs so that we don't need to apply the change of basis for the exact dirichlet trick in PreSolve */ 1776 if (pcbddc->benign_apply_coarse_only && pcbddc->use_exact_dirichlet_trick && pcbddc->change_interior) { 1777 PetscCall(VecCopy(r, pcis->vec1_global)); 1778 PetscCall(VecLockReadPush(pcis->vec1_global)); 1779 } 1780 } 1781 if (pcbddc->benign_have_null) { /* get p0 from r */ 1782 PetscCall(PCBDDCBenignGetOrSetP0(pc, r, PETSC_TRUE)); 1783 } 1784 if (pcbddc->interface_extension == PC_BDDC_INTERFACE_EXT_DIRICHLET && !pcbddc->exact_dirichlet_trick_app && !pcbddc->benign_apply_coarse_only) { 1785 PetscCall(VecCopy(r, z)); 1786 /* First Dirichlet solve */ 1787 PetscCall(VecScatterBegin(pcis->global_to_D, r, pcis->vec1_D, INSERT_VALUES, SCATTER_FORWARD)); 1788 PetscCall(VecScatterEnd(pcis->global_to_D, r, pcis->vec1_D, INSERT_VALUES, SCATTER_FORWARD)); 1789 /* 1790 Assembling right hand side for BDDC operator 1791 - pcis->vec1_D for the Dirichlet part (if needed, i.e. pcbddc->switch_static == PETSC_TRUE) 1792 - pcis->vec1_B the interface part of the global vector z 1793 */ 1794 if (n_D) { 1795 PetscCall(PetscLogEventBegin(PC_BDDC_Solves[pcbddc->current_level][0], pc, 0, 0, 0)); 1796 PetscCall(KSPSolve(pcbddc->ksp_D, pcis->vec1_D, pcis->vec2_D)); 1797 PetscCall(PetscLogEventEnd(PC_BDDC_Solves[pcbddc->current_level][0], pc, 0, 0, 0)); 1798 PetscCall(KSPCheckSolve(pcbddc->ksp_D, pc, pcis->vec2_D)); 1799 PetscCall(VecScale(pcis->vec2_D, m_one)); 1800 if (pcbddc->switch_static) { 1801 PetscCall(VecSet(pcis->vec1_N, 0.)); 1802 PetscCall(VecScatterBegin(pcis->N_to_D, pcis->vec2_D, pcis->vec1_N, INSERT_VALUES, SCATTER_REVERSE)); 1803 PetscCall(VecScatterEnd(pcis->N_to_D, pcis->vec2_D, pcis->vec1_N, INSERT_VALUES, SCATTER_REVERSE)); 1804 if (!pcbddc->switch_static_change) { 1805 PetscCall(MatMult(lA, pcis->vec1_N, pcis->vec2_N)); 1806 } else { 1807 PetscCall(MatMult(pcbddc->switch_static_change, pcis->vec1_N, pcis->vec2_N)); 1808 PetscCall(MatMult(lA, pcis->vec2_N, pcis->vec1_N)); 1809 PetscCall(MatMultTranspose(pcbddc->switch_static_change, pcis->vec1_N, pcis->vec2_N)); 1810 } 1811 PetscCall(VecScatterBegin(pcis->N_to_D, pcis->vec2_N, pcis->vec1_D, ADD_VALUES, SCATTER_FORWARD)); 1812 PetscCall(VecScatterEnd(pcis->N_to_D, pcis->vec2_N, pcis->vec1_D, ADD_VALUES, SCATTER_FORWARD)); 1813 PetscCall(VecScatterBegin(pcis->N_to_B, pcis->vec2_N, pcis->vec1_B, INSERT_VALUES, SCATTER_FORWARD)); 1814 PetscCall(VecScatterEnd(pcis->N_to_B, pcis->vec2_N, pcis->vec1_B, INSERT_VALUES, SCATTER_FORWARD)); 1815 } else { 1816 PetscCall(MatMult(pcis->A_BI, pcis->vec2_D, pcis->vec1_B)); 1817 } 1818 } else { 1819 PetscCall(VecSet(pcis->vec1_B, zero)); 1820 } 1821 PetscCall(VecScatterBegin(pcis->global_to_B, pcis->vec1_B, z, ADD_VALUES, SCATTER_REVERSE)); 1822 PetscCall(VecScatterEnd(pcis->global_to_B, pcis->vec1_B, z, ADD_VALUES, SCATTER_REVERSE)); 1823 PetscCall(PCBDDCScalingRestriction(pc, z, pcis->vec1_B)); 1824 } else { 1825 if (!pcbddc->benign_apply_coarse_only) PetscCall(PCBDDCScalingRestriction(pc, r, pcis->vec1_B)); 1826 } 1827 if (pcbddc->interface_extension == PC_BDDC_INTERFACE_EXT_LUMP) { 1828 PetscCheck(pcbddc->switch_static, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "You forgot to pass -pc_bddc_switch_static"); 1829 PetscCall(VecScatterBegin(pcis->global_to_D, r, pcis->vec1_D, INSERT_VALUES, SCATTER_FORWARD)); 1830 PetscCall(VecScatterEnd(pcis->global_to_D, r, pcis->vec1_D, INSERT_VALUES, SCATTER_FORWARD)); 1831 } 1832 1833 /* Apply interface preconditioner 1834 input/output vecs: pcis->vec1_B and pcis->vec1_D */ 1835 PetscCall(PCBDDCApplyInterfacePreconditioner(pc, PETSC_FALSE)); 1836 1837 /* Apply transpose of partition of unity operator */ 1838 PetscCall(PCBDDCScalingExtension(pc, pcis->vec1_B, z)); 1839 if (pcbddc->interface_extension == PC_BDDC_INTERFACE_EXT_LUMP) { 1840 PetscCall(VecScatterBegin(pcis->global_to_D, pcis->vec1_D, z, INSERT_VALUES, SCATTER_REVERSE)); 1841 PetscCall(VecScatterEnd(pcis->global_to_D, pcis->vec1_D, z, INSERT_VALUES, SCATTER_REVERSE)); 1842 PetscFunctionReturn(0); 1843 } 1844 /* Second Dirichlet solve and assembling of output */ 1845 PetscCall(VecScatterBegin(pcis->global_to_B, z, pcis->vec1_B, INSERT_VALUES, SCATTER_FORWARD)); 1846 PetscCall(VecScatterEnd(pcis->global_to_B, z, pcis->vec1_B, INSERT_VALUES, SCATTER_FORWARD)); 1847 if (n_B) { 1848 if (pcbddc->switch_static) { 1849 PetscCall(VecScatterBegin(pcis->N_to_D, pcis->vec1_D, pcis->vec1_N, INSERT_VALUES, SCATTER_REVERSE)); 1850 PetscCall(VecScatterEnd(pcis->N_to_D, pcis->vec1_D, pcis->vec1_N, INSERT_VALUES, SCATTER_REVERSE)); 1851 PetscCall(VecScatterBegin(pcis->N_to_B, pcis->vec1_B, pcis->vec1_N, INSERT_VALUES, SCATTER_REVERSE)); 1852 PetscCall(VecScatterEnd(pcis->N_to_B, pcis->vec1_B, pcis->vec1_N, INSERT_VALUES, SCATTER_REVERSE)); 1853 if (!pcbddc->switch_static_change) { 1854 PetscCall(MatMult(lA, pcis->vec1_N, pcis->vec2_N)); 1855 } else { 1856 PetscCall(MatMult(pcbddc->switch_static_change, pcis->vec1_N, pcis->vec2_N)); 1857 PetscCall(MatMult(lA, pcis->vec2_N, pcis->vec1_N)); 1858 PetscCall(MatMultTranspose(pcbddc->switch_static_change, pcis->vec1_N, pcis->vec2_N)); 1859 } 1860 PetscCall(VecScatterBegin(pcis->N_to_D, pcis->vec2_N, pcis->vec3_D, INSERT_VALUES, SCATTER_FORWARD)); 1861 PetscCall(VecScatterEnd(pcis->N_to_D, pcis->vec2_N, pcis->vec3_D, INSERT_VALUES, SCATTER_FORWARD)); 1862 } else { 1863 PetscCall(MatMult(pcis->A_IB, pcis->vec1_B, pcis->vec3_D)); 1864 } 1865 } else if (pcbddc->switch_static) { /* n_B is zero */ 1866 if (!pcbddc->switch_static_change) { 1867 PetscCall(MatMult(lA, pcis->vec1_D, pcis->vec3_D)); 1868 } else { 1869 PetscCall(MatMult(pcbddc->switch_static_change, pcis->vec1_D, pcis->vec1_N)); 1870 PetscCall(MatMult(lA, pcis->vec1_N, pcis->vec2_N)); 1871 PetscCall(MatMultTranspose(pcbddc->switch_static_change, pcis->vec2_N, pcis->vec3_D)); 1872 } 1873 } 1874 PetscCall(PetscLogEventBegin(PC_BDDC_Solves[pcbddc->current_level][0], pc, 0, 0, 0)); 1875 PetscCall(KSPSolve(pcbddc->ksp_D, pcis->vec3_D, pcis->vec4_D)); 1876 PetscCall(PetscLogEventEnd(PC_BDDC_Solves[pcbddc->current_level][0], pc, 0, 0, 0)); 1877 PetscCall(KSPCheckSolve(pcbddc->ksp_D, pc, pcis->vec4_D)); 1878 1879 if (!pcbddc->exact_dirichlet_trick_app && !pcbddc->benign_apply_coarse_only) { 1880 if (pcbddc->switch_static) { 1881 PetscCall(VecAXPBYPCZ(pcis->vec2_D, m_one, one, m_one, pcis->vec4_D, pcis->vec1_D)); 1882 } else { 1883 PetscCall(VecAXPBY(pcis->vec2_D, m_one, m_one, pcis->vec4_D)); 1884 } 1885 PetscCall(VecScatterBegin(pcis->global_to_D, pcis->vec2_D, z, INSERT_VALUES, SCATTER_REVERSE)); 1886 PetscCall(VecScatterEnd(pcis->global_to_D, pcis->vec2_D, z, INSERT_VALUES, SCATTER_REVERSE)); 1887 } else { 1888 if (pcbddc->switch_static) { 1889 PetscCall(VecAXPBY(pcis->vec4_D, one, m_one, pcis->vec1_D)); 1890 } else { 1891 PetscCall(VecScale(pcis->vec4_D, m_one)); 1892 } 1893 PetscCall(VecScatterBegin(pcis->global_to_D, pcis->vec4_D, z, INSERT_VALUES, SCATTER_REVERSE)); 1894 PetscCall(VecScatterEnd(pcis->global_to_D, pcis->vec4_D, z, INSERT_VALUES, SCATTER_REVERSE)); 1895 } 1896 if (pcbddc->benign_have_null) { /* set p0 (computed in PCBDDCApplyInterface) */ 1897 if (pcbddc->benign_apply_coarse_only) PetscCall(PetscArrayzero(pcbddc->benign_p0, pcbddc->benign_n)); 1898 PetscCall(PCBDDCBenignGetOrSetP0(pc, z, PETSC_FALSE)); 1899 } 1900 1901 if (pcbddc->ChangeOfBasisMatrix) { 1902 pcbddc->work_change = r; 1903 PetscCall(VecCopy(z, pcbddc->work_change)); 1904 PetscCall(MatMult(pcbddc->ChangeOfBasisMatrix, pcbddc->work_change, z)); 1905 } 1906 PetscFunctionReturn(0); 1907 } 1908 1909 /* 1910 PCApplyTranspose_BDDC - Applies the transpose of the BDDC operator to a vector. 1911 1912 Input Parameters: 1913 + pc - the preconditioner context 1914 - r - input vector (global) 1915 1916 Output Parameter: 1917 . z - output vector (global) 1918 1919 Application Interface Routine: PCApplyTranspose() 1920 */ 1921 PetscErrorCode PCApplyTranspose_BDDC(PC pc, Vec r, Vec z) 1922 { 1923 PC_IS *pcis = (PC_IS *)(pc->data); 1924 PC_BDDC *pcbddc = (PC_BDDC *)(pc->data); 1925 Mat lA = NULL; 1926 PetscInt n_B = pcis->n_B, n_D = pcis->n - n_B; 1927 const PetscScalar one = 1.0; 1928 const PetscScalar m_one = -1.0; 1929 const PetscScalar zero = 0.0; 1930 1931 PetscFunctionBegin; 1932 PetscCall(PetscCitationsRegister(citation, &cited)); 1933 if (pcbddc->switch_static) PetscCall(MatISGetLocalMat(pc->useAmat ? pc->mat : pc->pmat, &lA)); 1934 if (pcbddc->ChangeOfBasisMatrix) { 1935 Vec swap; 1936 1937 PetscCall(MatMultTranspose(pcbddc->ChangeOfBasisMatrix, r, pcbddc->work_change)); 1938 swap = pcbddc->work_change; 1939 pcbddc->work_change = r; 1940 r = swap; 1941 /* save rhs so that we don't need to apply the change of basis for the exact dirichlet trick in PreSolve */ 1942 if (pcbddc->benign_apply_coarse_only && pcbddc->exact_dirichlet_trick_app && pcbddc->change_interior) { 1943 PetscCall(VecCopy(r, pcis->vec1_global)); 1944 PetscCall(VecLockReadPush(pcis->vec1_global)); 1945 } 1946 } 1947 if (pcbddc->benign_have_null) { /* get p0 from r */ 1948 PetscCall(PCBDDCBenignGetOrSetP0(pc, r, PETSC_TRUE)); 1949 } 1950 if (!pcbddc->exact_dirichlet_trick_app && !pcbddc->benign_apply_coarse_only) { 1951 PetscCall(VecCopy(r, z)); 1952 /* First Dirichlet solve */ 1953 PetscCall(VecScatterBegin(pcis->global_to_D, r, pcis->vec1_D, INSERT_VALUES, SCATTER_FORWARD)); 1954 PetscCall(VecScatterEnd(pcis->global_to_D, r, pcis->vec1_D, INSERT_VALUES, SCATTER_FORWARD)); 1955 /* 1956 Assembling right hand side for BDDC operator 1957 - pcis->vec1_D for the Dirichlet part (if needed, i.e. pcbddc->switch_static == PETSC_TRUE) 1958 - pcis->vec1_B the interface part of the global vector z 1959 */ 1960 if (n_D) { 1961 PetscCall(PetscLogEventBegin(PC_BDDC_Solves[pcbddc->current_level][0], pc, 0, 0, 0)); 1962 PetscCall(KSPSolveTranspose(pcbddc->ksp_D, pcis->vec1_D, pcis->vec2_D)); 1963 PetscCall(PetscLogEventEnd(PC_BDDC_Solves[pcbddc->current_level][0], pc, 0, 0, 0)); 1964 PetscCall(KSPCheckSolve(pcbddc->ksp_D, pc, pcis->vec2_D)); 1965 PetscCall(VecScale(pcis->vec2_D, m_one)); 1966 if (pcbddc->switch_static) { 1967 PetscCall(VecSet(pcis->vec1_N, 0.)); 1968 PetscCall(VecScatterBegin(pcis->N_to_D, pcis->vec2_D, pcis->vec1_N, INSERT_VALUES, SCATTER_REVERSE)); 1969 PetscCall(VecScatterEnd(pcis->N_to_D, pcis->vec2_D, pcis->vec1_N, INSERT_VALUES, SCATTER_REVERSE)); 1970 if (!pcbddc->switch_static_change) { 1971 PetscCall(MatMultTranspose(lA, pcis->vec1_N, pcis->vec2_N)); 1972 } else { 1973 PetscCall(MatMult(pcbddc->switch_static_change, pcis->vec1_N, pcis->vec2_N)); 1974 PetscCall(MatMultTranspose(lA, pcis->vec2_N, pcis->vec1_N)); 1975 PetscCall(MatMultTranspose(pcbddc->switch_static_change, pcis->vec1_N, pcis->vec2_N)); 1976 } 1977 PetscCall(VecScatterBegin(pcis->N_to_D, pcis->vec2_N, pcis->vec1_D, ADD_VALUES, SCATTER_FORWARD)); 1978 PetscCall(VecScatterEnd(pcis->N_to_D, pcis->vec2_N, pcis->vec1_D, ADD_VALUES, SCATTER_FORWARD)); 1979 PetscCall(VecScatterBegin(pcis->N_to_B, pcis->vec2_N, pcis->vec1_B, INSERT_VALUES, SCATTER_FORWARD)); 1980 PetscCall(VecScatterEnd(pcis->N_to_B, pcis->vec2_N, pcis->vec1_B, INSERT_VALUES, SCATTER_FORWARD)); 1981 } else { 1982 PetscCall(MatMultTranspose(pcis->A_IB, pcis->vec2_D, pcis->vec1_B)); 1983 } 1984 } else { 1985 PetscCall(VecSet(pcis->vec1_B, zero)); 1986 } 1987 PetscCall(VecScatterBegin(pcis->global_to_B, pcis->vec1_B, z, ADD_VALUES, SCATTER_REVERSE)); 1988 PetscCall(VecScatterEnd(pcis->global_to_B, pcis->vec1_B, z, ADD_VALUES, SCATTER_REVERSE)); 1989 PetscCall(PCBDDCScalingRestriction(pc, z, pcis->vec1_B)); 1990 } else { 1991 PetscCall(PCBDDCScalingRestriction(pc, r, pcis->vec1_B)); 1992 } 1993 1994 /* Apply interface preconditioner 1995 input/output vecs: pcis->vec1_B and pcis->vec1_D */ 1996 PetscCall(PCBDDCApplyInterfacePreconditioner(pc, PETSC_TRUE)); 1997 1998 /* Apply transpose of partition of unity operator */ 1999 PetscCall(PCBDDCScalingExtension(pc, pcis->vec1_B, z)); 2000 2001 /* Second Dirichlet solve and assembling of output */ 2002 PetscCall(VecScatterBegin(pcis->global_to_B, z, pcis->vec1_B, INSERT_VALUES, SCATTER_FORWARD)); 2003 PetscCall(VecScatterEnd(pcis->global_to_B, z, pcis->vec1_B, INSERT_VALUES, SCATTER_FORWARD)); 2004 if (n_B) { 2005 if (pcbddc->switch_static) { 2006 PetscCall(VecScatterBegin(pcis->N_to_D, pcis->vec1_D, pcis->vec1_N, INSERT_VALUES, SCATTER_REVERSE)); 2007 PetscCall(VecScatterEnd(pcis->N_to_D, pcis->vec1_D, pcis->vec1_N, INSERT_VALUES, SCATTER_REVERSE)); 2008 PetscCall(VecScatterBegin(pcis->N_to_B, pcis->vec1_B, pcis->vec1_N, INSERT_VALUES, SCATTER_REVERSE)); 2009 PetscCall(VecScatterEnd(pcis->N_to_B, pcis->vec1_B, pcis->vec1_N, INSERT_VALUES, SCATTER_REVERSE)); 2010 if (!pcbddc->switch_static_change) { 2011 PetscCall(MatMultTranspose(lA, pcis->vec1_N, pcis->vec2_N)); 2012 } else { 2013 PetscCall(MatMult(pcbddc->switch_static_change, pcis->vec1_N, pcis->vec2_N)); 2014 PetscCall(MatMultTranspose(lA, pcis->vec2_N, pcis->vec1_N)); 2015 PetscCall(MatMultTranspose(pcbddc->switch_static_change, pcis->vec1_N, pcis->vec2_N)); 2016 } 2017 PetscCall(VecScatterBegin(pcis->N_to_D, pcis->vec2_N, pcis->vec3_D, INSERT_VALUES, SCATTER_FORWARD)); 2018 PetscCall(VecScatterEnd(pcis->N_to_D, pcis->vec2_N, pcis->vec3_D, INSERT_VALUES, SCATTER_FORWARD)); 2019 } else { 2020 PetscCall(MatMultTranspose(pcis->A_BI, pcis->vec1_B, pcis->vec3_D)); 2021 } 2022 } else if (pcbddc->switch_static) { /* n_B is zero */ 2023 if (!pcbddc->switch_static_change) { 2024 PetscCall(MatMultTranspose(lA, pcis->vec1_D, pcis->vec3_D)); 2025 } else { 2026 PetscCall(MatMult(pcbddc->switch_static_change, pcis->vec1_D, pcis->vec1_N)); 2027 PetscCall(MatMultTranspose(lA, pcis->vec1_N, pcis->vec2_N)); 2028 PetscCall(MatMultTranspose(pcbddc->switch_static_change, pcis->vec2_N, pcis->vec3_D)); 2029 } 2030 } 2031 PetscCall(PetscLogEventBegin(PC_BDDC_Solves[pcbddc->current_level][0], pc, 0, 0, 0)); 2032 PetscCall(KSPSolveTranspose(pcbddc->ksp_D, pcis->vec3_D, pcis->vec4_D)); 2033 PetscCall(PetscLogEventEnd(PC_BDDC_Solves[pcbddc->current_level][0], pc, 0, 0, 0)); 2034 PetscCall(KSPCheckSolve(pcbddc->ksp_D, pc, pcis->vec4_D)); 2035 if (!pcbddc->exact_dirichlet_trick_app && !pcbddc->benign_apply_coarse_only) { 2036 if (pcbddc->switch_static) { 2037 PetscCall(VecAXPBYPCZ(pcis->vec2_D, m_one, one, m_one, pcis->vec4_D, pcis->vec1_D)); 2038 } else { 2039 PetscCall(VecAXPBY(pcis->vec2_D, m_one, m_one, pcis->vec4_D)); 2040 } 2041 PetscCall(VecScatterBegin(pcis->global_to_D, pcis->vec2_D, z, INSERT_VALUES, SCATTER_REVERSE)); 2042 PetscCall(VecScatterEnd(pcis->global_to_D, pcis->vec2_D, z, INSERT_VALUES, SCATTER_REVERSE)); 2043 } else { 2044 if (pcbddc->switch_static) { 2045 PetscCall(VecAXPBY(pcis->vec4_D, one, m_one, pcis->vec1_D)); 2046 } else { 2047 PetscCall(VecScale(pcis->vec4_D, m_one)); 2048 } 2049 PetscCall(VecScatterBegin(pcis->global_to_D, pcis->vec4_D, z, INSERT_VALUES, SCATTER_REVERSE)); 2050 PetscCall(VecScatterEnd(pcis->global_to_D, pcis->vec4_D, z, INSERT_VALUES, SCATTER_REVERSE)); 2051 } 2052 if (pcbddc->benign_have_null) { /* set p0 (computed in PCBDDCApplyInterface) */ 2053 PetscCall(PCBDDCBenignGetOrSetP0(pc, z, PETSC_FALSE)); 2054 } 2055 if (pcbddc->ChangeOfBasisMatrix) { 2056 pcbddc->work_change = r; 2057 PetscCall(VecCopy(z, pcbddc->work_change)); 2058 PetscCall(MatMult(pcbddc->ChangeOfBasisMatrix, pcbddc->work_change, z)); 2059 } 2060 PetscFunctionReturn(0); 2061 } 2062 2063 PetscErrorCode PCReset_BDDC(PC pc) 2064 { 2065 PC_BDDC *pcbddc = (PC_BDDC *)pc->data; 2066 PC_IS *pcis = (PC_IS *)pc->data; 2067 KSP kspD, kspR, kspC; 2068 2069 PetscFunctionBegin; 2070 /* free BDDC custom data */ 2071 PetscCall(PCBDDCResetCustomization(pc)); 2072 /* destroy objects related to topography */ 2073 PetscCall(PCBDDCResetTopography(pc)); 2074 /* destroy objects for scaling operator */ 2075 PetscCall(PCBDDCScalingDestroy(pc)); 2076 /* free solvers stuff */ 2077 PetscCall(PCBDDCResetSolvers(pc)); 2078 /* free global vectors needed in presolve */ 2079 PetscCall(VecDestroy(&pcbddc->temp_solution)); 2080 PetscCall(VecDestroy(&pcbddc->original_rhs)); 2081 /* free data created by PCIS */ 2082 PetscCall(PCISDestroy(pc)); 2083 2084 /* restore defaults */ 2085 kspD = pcbddc->ksp_D; 2086 kspR = pcbddc->ksp_R; 2087 kspC = pcbddc->coarse_ksp; 2088 PetscCall(PetscMemzero(pc->data, sizeof(*pcbddc))); 2089 pcis->n_neigh = -1; 2090 pcis->scaling_factor = 1.0; 2091 pcis->reusesubmatrices = PETSC_TRUE; 2092 pcbddc->use_local_adj = PETSC_TRUE; 2093 pcbddc->use_vertices = PETSC_TRUE; 2094 pcbddc->use_edges = PETSC_TRUE; 2095 pcbddc->symmetric_primal = PETSC_TRUE; 2096 pcbddc->vertex_size = 1; 2097 pcbddc->recompute_topography = PETSC_TRUE; 2098 pcbddc->coarse_size = -1; 2099 pcbddc->use_exact_dirichlet_trick = PETSC_TRUE; 2100 pcbddc->coarsening_ratio = 8; 2101 pcbddc->coarse_eqs_per_proc = 1; 2102 pcbddc->benign_compute_correction = PETSC_TRUE; 2103 pcbddc->nedfield = -1; 2104 pcbddc->nedglobal = PETSC_TRUE; 2105 pcbddc->graphmaxcount = PETSC_MAX_INT; 2106 pcbddc->sub_schurs_layers = -1; 2107 pcbddc->ksp_D = kspD; 2108 pcbddc->ksp_R = kspR; 2109 pcbddc->coarse_ksp = kspC; 2110 PetscFunctionReturn(0); 2111 } 2112 2113 PetscErrorCode PCDestroy_BDDC(PC pc) 2114 { 2115 PC_BDDC *pcbddc = (PC_BDDC *)pc->data; 2116 2117 PetscFunctionBegin; 2118 PetscCall(PCReset_BDDC(pc)); 2119 PetscCall(KSPDestroy(&pcbddc->ksp_D)); 2120 PetscCall(KSPDestroy(&pcbddc->ksp_R)); 2121 PetscCall(KSPDestroy(&pcbddc->coarse_ksp)); 2122 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBDDCSetDiscreteGradient_C", NULL)); 2123 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBDDCSetDivergenceMat_C", NULL)); 2124 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBDDCSetChangeOfBasisMat_C", NULL)); 2125 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBDDCSetPrimalVerticesLocalIS_C", NULL)); 2126 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBDDCSetPrimalVerticesIS_C", NULL)); 2127 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBDDCGetPrimalVerticesLocalIS_C", NULL)); 2128 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBDDCGetPrimalVerticesIS_C", NULL)); 2129 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBDDCSetCoarseningRatio_C", NULL)); 2130 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBDDCSetLevel_C", NULL)); 2131 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBDDCSetUseExactDirichlet_C", NULL)); 2132 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBDDCSetLevels_C", NULL)); 2133 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBDDCSetDirichletBoundaries_C", NULL)); 2134 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBDDCSetDirichletBoundariesLocal_C", NULL)); 2135 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBDDCSetNeumannBoundaries_C", NULL)); 2136 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBDDCSetNeumannBoundariesLocal_C", NULL)); 2137 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBDDCGetDirichletBoundaries_C", NULL)); 2138 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBDDCGetDirichletBoundariesLocal_C", NULL)); 2139 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBDDCGetNeumannBoundaries_C", NULL)); 2140 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBDDCGetNeumannBoundariesLocal_C", NULL)); 2141 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBDDCSetDofsSplitting_C", NULL)); 2142 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBDDCSetDofsSplittingLocal_C", NULL)); 2143 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBDDCSetLocalAdjacencyGraph_C", NULL)); 2144 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBDDCCreateFETIDPOperators_C", NULL)); 2145 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBDDCMatFETIDPGetRHS_C", NULL)); 2146 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBDDCMatFETIDPGetSolution_C", NULL)); 2147 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCPreSolveChangeRHS_C", NULL)); 2148 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSetCoordinates_C", NULL)); 2149 PetscCall(PetscFree(pc->data)); 2150 PetscFunctionReturn(0); 2151 } 2152 2153 static PetscErrorCode PCSetCoordinates_BDDC(PC pc, PetscInt dim, PetscInt nloc, PetscReal *coords) 2154 { 2155 PC_BDDC *pcbddc = (PC_BDDC *)pc->data; 2156 PCBDDCGraph mat_graph = pcbddc->mat_graph; 2157 2158 PetscFunctionBegin; 2159 PetscCall(PetscFree(mat_graph->coords)); 2160 PetscCall(PetscMalloc1(nloc * dim, &mat_graph->coords)); 2161 PetscCall(PetscArraycpy(mat_graph->coords, coords, nloc * dim)); 2162 mat_graph->cnloc = nloc; 2163 mat_graph->cdim = dim; 2164 mat_graph->cloc = PETSC_FALSE; 2165 /* flg setup */ 2166 pcbddc->recompute_topography = PETSC_TRUE; 2167 pcbddc->corner_selected = PETSC_FALSE; 2168 PetscFunctionReturn(0); 2169 } 2170 2171 static PetscErrorCode PCPreSolveChangeRHS_BDDC(PC pc, PetscBool *change) 2172 { 2173 PetscFunctionBegin; 2174 *change = PETSC_TRUE; 2175 PetscFunctionReturn(0); 2176 } 2177 2178 static PetscErrorCode PCBDDCMatFETIDPGetRHS_BDDC(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs) 2179 { 2180 FETIDPMat_ctx mat_ctx; 2181 Vec work; 2182 PC_IS *pcis; 2183 PC_BDDC *pcbddc; 2184 2185 PetscFunctionBegin; 2186 PetscCall(MatShellGetContext(fetidp_mat, &mat_ctx)); 2187 pcis = (PC_IS *)mat_ctx->pc->data; 2188 pcbddc = (PC_BDDC *)mat_ctx->pc->data; 2189 2190 PetscCall(VecSet(fetidp_flux_rhs, 0.0)); 2191 /* copy rhs since we may change it during PCPreSolve_BDDC */ 2192 if (!pcbddc->original_rhs) PetscCall(VecDuplicate(pcis->vec1_global, &pcbddc->original_rhs)); 2193 if (mat_ctx->rhs_flip) { 2194 PetscCall(VecPointwiseMult(pcbddc->original_rhs, standard_rhs, mat_ctx->rhs_flip)); 2195 } else { 2196 PetscCall(VecCopy(standard_rhs, pcbddc->original_rhs)); 2197 } 2198 if (mat_ctx->g2g_p) { 2199 /* interface pressure rhs */ 2200 PetscCall(VecScatterBegin(mat_ctx->g2g_p, fetidp_flux_rhs, pcbddc->original_rhs, INSERT_VALUES, SCATTER_REVERSE)); 2201 PetscCall(VecScatterEnd(mat_ctx->g2g_p, fetidp_flux_rhs, pcbddc->original_rhs, INSERT_VALUES, SCATTER_REVERSE)); 2202 PetscCall(VecScatterBegin(mat_ctx->g2g_p, standard_rhs, fetidp_flux_rhs, INSERT_VALUES, SCATTER_FORWARD)); 2203 PetscCall(VecScatterEnd(mat_ctx->g2g_p, standard_rhs, fetidp_flux_rhs, INSERT_VALUES, SCATTER_FORWARD)); 2204 if (!mat_ctx->rhs_flip) PetscCall(VecScale(fetidp_flux_rhs, -1.)); 2205 } 2206 /* 2207 change of basis for physical rhs if needed 2208 It also changes the rhs in case of dirichlet boundaries 2209 */ 2210 PetscCall(PCPreSolve_BDDC(mat_ctx->pc, NULL, pcbddc->original_rhs, NULL)); 2211 if (pcbddc->ChangeOfBasisMatrix) { 2212 PetscCall(MatMultTranspose(pcbddc->ChangeOfBasisMatrix, pcbddc->original_rhs, pcbddc->work_change)); 2213 work = pcbddc->work_change; 2214 } else { 2215 work = pcbddc->original_rhs; 2216 } 2217 /* store vectors for computation of fetidp final solution */ 2218 PetscCall(VecScatterBegin(pcis->global_to_D, work, mat_ctx->temp_solution_D, INSERT_VALUES, SCATTER_FORWARD)); 2219 PetscCall(VecScatterEnd(pcis->global_to_D, work, mat_ctx->temp_solution_D, INSERT_VALUES, SCATTER_FORWARD)); 2220 /* scale rhs since it should be unassembled */ 2221 /* TODO use counter scaling? (also below) */ 2222 PetscCall(VecScatterBegin(pcis->global_to_B, work, mat_ctx->temp_solution_B, INSERT_VALUES, SCATTER_FORWARD)); 2223 PetscCall(VecScatterEnd(pcis->global_to_B, work, mat_ctx->temp_solution_B, INSERT_VALUES, SCATTER_FORWARD)); 2224 /* Apply partition of unity */ 2225 PetscCall(VecPointwiseMult(mat_ctx->temp_solution_B, pcis->D, mat_ctx->temp_solution_B)); 2226 /* PetscCall(PCBDDCScalingRestriction(mat_ctx->pc,work,mat_ctx->temp_solution_B)); */ 2227 if (!pcbddc->switch_static) { 2228 /* compute partially subassembled Schur complement right-hand side */ 2229 PetscCall(PetscLogEventBegin(PC_BDDC_Solves[pcbddc->current_level][0], mat_ctx->pc, 0, 0, 0)); 2230 PetscCall(KSPSolve(pcbddc->ksp_D, mat_ctx->temp_solution_D, pcis->vec1_D)); 2231 PetscCall(PetscLogEventEnd(PC_BDDC_Solves[pcbddc->current_level][0], mat_ctx->pc, 0, 0, 0)); 2232 /* Cannot propagate up error in KSPSolve() because there is no access to the PC */ 2233 PetscCall(MatMult(pcis->A_BI, pcis->vec1_D, pcis->vec1_B)); 2234 PetscCall(VecAXPY(mat_ctx->temp_solution_B, -1.0, pcis->vec1_B)); 2235 PetscCall(VecSet(work, 0.0)); 2236 PetscCall(VecScatterBegin(pcis->global_to_B, mat_ctx->temp_solution_B, work, ADD_VALUES, SCATTER_REVERSE)); 2237 PetscCall(VecScatterEnd(pcis->global_to_B, mat_ctx->temp_solution_B, work, ADD_VALUES, SCATTER_REVERSE)); 2238 /* PetscCall(PCBDDCScalingRestriction(mat_ctx->pc,work,mat_ctx->temp_solution_B)); */ 2239 PetscCall(VecScatterBegin(pcis->global_to_B, work, mat_ctx->temp_solution_B, INSERT_VALUES, SCATTER_FORWARD)); 2240 PetscCall(VecScatterEnd(pcis->global_to_B, work, mat_ctx->temp_solution_B, INSERT_VALUES, SCATTER_FORWARD)); 2241 PetscCall(VecPointwiseMult(mat_ctx->temp_solution_B, pcis->D, mat_ctx->temp_solution_B)); 2242 } 2243 /* BDDC rhs */ 2244 PetscCall(VecCopy(mat_ctx->temp_solution_B, pcis->vec1_B)); 2245 if (pcbddc->switch_static) PetscCall(VecCopy(mat_ctx->temp_solution_D, pcis->vec1_D)); 2246 /* apply BDDC */ 2247 PetscCall(PetscArrayzero(pcbddc->benign_p0, pcbddc->benign_n)); 2248 PetscCall(PCBDDCApplyInterfacePreconditioner(mat_ctx->pc, PETSC_FALSE)); 2249 PetscCall(PetscArrayzero(pcbddc->benign_p0, pcbddc->benign_n)); 2250 2251 /* Application of B_delta and assembling of rhs for fetidp fluxes */ 2252 PetscCall(MatMult(mat_ctx->B_delta, pcis->vec1_B, mat_ctx->lambda_local)); 2253 PetscCall(VecScatterBegin(mat_ctx->l2g_lambda, mat_ctx->lambda_local, fetidp_flux_rhs, ADD_VALUES, SCATTER_FORWARD)); 2254 PetscCall(VecScatterEnd(mat_ctx->l2g_lambda, mat_ctx->lambda_local, fetidp_flux_rhs, ADD_VALUES, SCATTER_FORWARD)); 2255 /* Add contribution to interface pressures */ 2256 if (mat_ctx->l2g_p) { 2257 PetscCall(MatMult(mat_ctx->B_BB, pcis->vec1_B, mat_ctx->vP)); 2258 if (pcbddc->switch_static) PetscCall(MatMultAdd(mat_ctx->B_BI, pcis->vec1_D, mat_ctx->vP, mat_ctx->vP)); 2259 PetscCall(VecScatterBegin(mat_ctx->l2g_p, mat_ctx->vP, fetidp_flux_rhs, ADD_VALUES, SCATTER_FORWARD)); 2260 PetscCall(VecScatterEnd(mat_ctx->l2g_p, mat_ctx->vP, fetidp_flux_rhs, ADD_VALUES, SCATTER_FORWARD)); 2261 } 2262 PetscFunctionReturn(0); 2263 } 2264 2265 /*@ 2266 PCBDDCMatFETIDPGetRHS - Compute the right-hand side for FETI-DP linear system using the physical right-hand side 2267 2268 Collective 2269 2270 Input Parameters: 2271 + fetidp_mat - the FETI-DP matrix object obtained by a call to `PCBDDCCreateFETIDPOperators()` 2272 - standard_rhs - the right-hand side of the original linear system 2273 2274 Output Parameter: 2275 . fetidp_flux_rhs - the right-hand side for the FETI-DP linear system 2276 2277 Level: developer 2278 2279 .seealso: `PCBDDC`, `PCBDDCCreateFETIDPOperators()`, `PCBDDCMatFETIDPGetSolution()` 2280 @*/ 2281 PetscErrorCode PCBDDCMatFETIDPGetRHS(Mat fetidp_mat, Vec standard_rhs, Vec fetidp_flux_rhs) 2282 { 2283 FETIDPMat_ctx mat_ctx; 2284 2285 PetscFunctionBegin; 2286 PetscValidHeaderSpecific(fetidp_mat, MAT_CLASSID, 1); 2287 PetscValidHeaderSpecific(standard_rhs, VEC_CLASSID, 2); 2288 PetscValidHeaderSpecific(fetidp_flux_rhs, VEC_CLASSID, 3); 2289 PetscCall(MatShellGetContext(fetidp_mat, &mat_ctx)); 2290 PetscUseMethod(mat_ctx->pc, "PCBDDCMatFETIDPGetRHS_C", (Mat, Vec, Vec), (fetidp_mat, standard_rhs, fetidp_flux_rhs)); 2291 PetscFunctionReturn(0); 2292 } 2293 2294 static PetscErrorCode PCBDDCMatFETIDPGetSolution_BDDC(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol) 2295 { 2296 FETIDPMat_ctx mat_ctx; 2297 PC_IS *pcis; 2298 PC_BDDC *pcbddc; 2299 Vec work; 2300 2301 PetscFunctionBegin; 2302 PetscCall(MatShellGetContext(fetidp_mat, &mat_ctx)); 2303 pcis = (PC_IS *)mat_ctx->pc->data; 2304 pcbddc = (PC_BDDC *)mat_ctx->pc->data; 2305 2306 /* apply B_delta^T */ 2307 PetscCall(VecSet(pcis->vec1_B, 0.)); 2308 PetscCall(VecScatterBegin(mat_ctx->l2g_lambda, fetidp_flux_sol, mat_ctx->lambda_local, INSERT_VALUES, SCATTER_REVERSE)); 2309 PetscCall(VecScatterEnd(mat_ctx->l2g_lambda, fetidp_flux_sol, mat_ctx->lambda_local, INSERT_VALUES, SCATTER_REVERSE)); 2310 PetscCall(MatMultTranspose(mat_ctx->B_delta, mat_ctx->lambda_local, pcis->vec1_B)); 2311 if (mat_ctx->l2g_p) { 2312 PetscCall(VecScatterBegin(mat_ctx->l2g_p, fetidp_flux_sol, mat_ctx->vP, INSERT_VALUES, SCATTER_REVERSE)); 2313 PetscCall(VecScatterEnd(mat_ctx->l2g_p, fetidp_flux_sol, mat_ctx->vP, INSERT_VALUES, SCATTER_REVERSE)); 2314 PetscCall(MatMultAdd(mat_ctx->Bt_BB, mat_ctx->vP, pcis->vec1_B, pcis->vec1_B)); 2315 } 2316 2317 /* compute rhs for BDDC application */ 2318 PetscCall(VecAYPX(pcis->vec1_B, -1.0, mat_ctx->temp_solution_B)); 2319 if (pcbddc->switch_static) { 2320 PetscCall(VecCopy(mat_ctx->temp_solution_D, pcis->vec1_D)); 2321 if (mat_ctx->l2g_p) { 2322 PetscCall(VecScale(mat_ctx->vP, -1.)); 2323 PetscCall(MatMultAdd(mat_ctx->Bt_BI, mat_ctx->vP, pcis->vec1_D, pcis->vec1_D)); 2324 } 2325 } 2326 2327 /* apply BDDC */ 2328 PetscCall(PetscArrayzero(pcbddc->benign_p0, pcbddc->benign_n)); 2329 PetscCall(PCBDDCApplyInterfacePreconditioner(mat_ctx->pc, PETSC_FALSE)); 2330 2331 /* put values into global vector */ 2332 if (pcbddc->ChangeOfBasisMatrix) work = pcbddc->work_change; 2333 else work = standard_sol; 2334 PetscCall(VecScatterBegin(pcis->global_to_B, pcis->vec1_B, work, INSERT_VALUES, SCATTER_REVERSE)); 2335 PetscCall(VecScatterEnd(pcis->global_to_B, pcis->vec1_B, work, INSERT_VALUES, SCATTER_REVERSE)); 2336 if (!pcbddc->switch_static) { 2337 /* compute values into the interior if solved for the partially subassembled Schur complement */ 2338 PetscCall(MatMult(pcis->A_IB, pcis->vec1_B, pcis->vec1_D)); 2339 PetscCall(VecAYPX(pcis->vec1_D, -1.0, mat_ctx->temp_solution_D)); 2340 PetscCall(PetscLogEventBegin(PC_BDDC_Solves[pcbddc->current_level][0], mat_ctx->pc, 0, 0, 0)); 2341 PetscCall(KSPSolve(pcbddc->ksp_D, pcis->vec1_D, pcis->vec1_D)); 2342 PetscCall(PetscLogEventEnd(PC_BDDC_Solves[pcbddc->current_level][0], mat_ctx->pc, 0, 0, 0)); 2343 /* Cannot propagate up error in KSPSolve() because there is no access to the PC */ 2344 } 2345 2346 PetscCall(VecScatterBegin(pcis->global_to_D, pcis->vec1_D, work, INSERT_VALUES, SCATTER_REVERSE)); 2347 PetscCall(VecScatterEnd(pcis->global_to_D, pcis->vec1_D, work, INSERT_VALUES, SCATTER_REVERSE)); 2348 /* add p0 solution to final solution */ 2349 PetscCall(PCBDDCBenignGetOrSetP0(mat_ctx->pc, work, PETSC_FALSE)); 2350 if (pcbddc->ChangeOfBasisMatrix) PetscCall(MatMult(pcbddc->ChangeOfBasisMatrix, work, standard_sol)); 2351 PetscCall(PCPostSolve_BDDC(mat_ctx->pc, NULL, NULL, standard_sol)); 2352 if (mat_ctx->g2g_p) { 2353 PetscCall(VecScatterBegin(mat_ctx->g2g_p, fetidp_flux_sol, standard_sol, INSERT_VALUES, SCATTER_REVERSE)); 2354 PetscCall(VecScatterEnd(mat_ctx->g2g_p, fetidp_flux_sol, standard_sol, INSERT_VALUES, SCATTER_REVERSE)); 2355 } 2356 PetscFunctionReturn(0); 2357 } 2358 2359 static PetscErrorCode PCView_BDDCIPC(PC pc, PetscViewer viewer) 2360 { 2361 BDDCIPC_ctx bddcipc_ctx; 2362 PetscBool isascii; 2363 2364 PetscFunctionBegin; 2365 PetscCall(PCShellGetContext(pc, &bddcipc_ctx)); 2366 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii)); 2367 if (isascii) PetscCall(PetscViewerASCIIPrintf(viewer, "BDDC interface preconditioner\n")); 2368 PetscCall(PetscViewerASCIIPushTab(viewer)); 2369 PetscCall(PCView(bddcipc_ctx->bddc, viewer)); 2370 PetscCall(PetscViewerASCIIPopTab(viewer)); 2371 PetscFunctionReturn(0); 2372 } 2373 2374 static PetscErrorCode PCSetUp_BDDCIPC(PC pc) 2375 { 2376 BDDCIPC_ctx bddcipc_ctx; 2377 PetscBool isbddc; 2378 Vec vv; 2379 IS is; 2380 PC_IS *pcis; 2381 2382 PetscFunctionBegin; 2383 PetscCall(PCShellGetContext(pc, &bddcipc_ctx)); 2384 PetscCall(PetscObjectTypeCompare((PetscObject)bddcipc_ctx->bddc, PCBDDC, &isbddc)); 2385 PetscCheck(isbddc, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Invalid type %s. Must be of type bddc", ((PetscObject)bddcipc_ctx->bddc)->type_name); 2386 PetscCall(PCSetUp(bddcipc_ctx->bddc)); 2387 2388 /* create interface scatter */ 2389 pcis = (PC_IS *)(bddcipc_ctx->bddc->data); 2390 PetscCall(VecScatterDestroy(&bddcipc_ctx->g2l)); 2391 PetscCall(MatCreateVecs(pc->pmat, &vv, NULL)); 2392 PetscCall(ISRenumber(pcis->is_B_global, NULL, NULL, &is)); 2393 PetscCall(VecScatterCreate(vv, is, pcis->vec1_B, NULL, &bddcipc_ctx->g2l)); 2394 PetscCall(ISDestroy(&is)); 2395 PetscCall(VecDestroy(&vv)); 2396 PetscFunctionReturn(0); 2397 } 2398 2399 static PetscErrorCode PCApply_BDDCIPC(PC pc, Vec r, Vec x) 2400 { 2401 BDDCIPC_ctx bddcipc_ctx; 2402 PC_IS *pcis; 2403 VecScatter tmps; 2404 2405 PetscFunctionBegin; 2406 PetscCall(PCShellGetContext(pc, &bddcipc_ctx)); 2407 pcis = (PC_IS *)(bddcipc_ctx->bddc->data); 2408 tmps = pcis->global_to_B; 2409 pcis->global_to_B = bddcipc_ctx->g2l; 2410 PetscCall(PCBDDCScalingRestriction(bddcipc_ctx->bddc, r, pcis->vec1_B)); 2411 PetscCall(PCBDDCApplyInterfacePreconditioner(bddcipc_ctx->bddc, PETSC_FALSE)); 2412 PetscCall(PCBDDCScalingExtension(bddcipc_ctx->bddc, pcis->vec1_B, x)); 2413 pcis->global_to_B = tmps; 2414 PetscFunctionReturn(0); 2415 } 2416 2417 static PetscErrorCode PCApplyTranspose_BDDCIPC(PC pc, Vec r, Vec x) 2418 { 2419 BDDCIPC_ctx bddcipc_ctx; 2420 PC_IS *pcis; 2421 VecScatter tmps; 2422 2423 PetscFunctionBegin; 2424 PetscCall(PCShellGetContext(pc, &bddcipc_ctx)); 2425 pcis = (PC_IS *)(bddcipc_ctx->bddc->data); 2426 tmps = pcis->global_to_B; 2427 pcis->global_to_B = bddcipc_ctx->g2l; 2428 PetscCall(PCBDDCScalingRestriction(bddcipc_ctx->bddc, r, pcis->vec1_B)); 2429 PetscCall(PCBDDCApplyInterfacePreconditioner(bddcipc_ctx->bddc, PETSC_TRUE)); 2430 PetscCall(PCBDDCScalingExtension(bddcipc_ctx->bddc, pcis->vec1_B, x)); 2431 pcis->global_to_B = tmps; 2432 PetscFunctionReturn(0); 2433 } 2434 2435 static PetscErrorCode PCDestroy_BDDCIPC(PC pc) 2436 { 2437 BDDCIPC_ctx bddcipc_ctx; 2438 2439 PetscFunctionBegin; 2440 PetscCall(PCShellGetContext(pc, &bddcipc_ctx)); 2441 PetscCall(PCDestroy(&bddcipc_ctx->bddc)); 2442 PetscCall(VecScatterDestroy(&bddcipc_ctx->g2l)); 2443 PetscCall(PetscFree(bddcipc_ctx)); 2444 PetscFunctionReturn(0); 2445 } 2446 2447 /*@ 2448 PCBDDCMatFETIDPGetSolution - Compute the physical solution using the solution of the FETI-DP linear system 2449 2450 Collective 2451 2452 Input Parameters: 2453 + fetidp_mat - the FETI-DP matrix obtained by a call to PCBDDCCreateFETIDPOperators()` 2454 - fetidp_flux_sol - the solution of the FETI-DP linear system` 2455 2456 Output Parameter: 2457 . standard_sol - the solution defined on the physical domain 2458 2459 Level: developer 2460 2461 .seealso: `PCBDDC`, `PCBDDCCreateFETIDPOperators()`, `PCBDDCMatFETIDPGetRHS()` 2462 @*/ 2463 PetscErrorCode PCBDDCMatFETIDPGetSolution(Mat fetidp_mat, Vec fetidp_flux_sol, Vec standard_sol) 2464 { 2465 FETIDPMat_ctx mat_ctx; 2466 2467 PetscFunctionBegin; 2468 PetscValidHeaderSpecific(fetidp_mat, MAT_CLASSID, 1); 2469 PetscValidHeaderSpecific(fetidp_flux_sol, VEC_CLASSID, 2); 2470 PetscValidHeaderSpecific(standard_sol, VEC_CLASSID, 3); 2471 PetscCall(MatShellGetContext(fetidp_mat, &mat_ctx)); 2472 PetscUseMethod(mat_ctx->pc, "PCBDDCMatFETIDPGetSolution_C", (Mat, Vec, Vec), (fetidp_mat, fetidp_flux_sol, standard_sol)); 2473 PetscFunctionReturn(0); 2474 } 2475 2476 static PetscErrorCode PCBDDCCreateFETIDPOperators_BDDC(PC pc, PetscBool fully_redundant, const char *prefix, Mat *fetidp_mat, PC *fetidp_pc) 2477 { 2478 FETIDPMat_ctx fetidpmat_ctx; 2479 Mat newmat; 2480 FETIDPPC_ctx fetidppc_ctx; 2481 PC newpc; 2482 MPI_Comm comm; 2483 2484 PetscFunctionBegin; 2485 PetscCall(PetscObjectGetComm((PetscObject)pc, &comm)); 2486 /* FETI-DP matrix */ 2487 PetscCall(PCBDDCCreateFETIDPMatContext(pc, &fetidpmat_ctx)); 2488 fetidpmat_ctx->fully_redundant = fully_redundant; 2489 PetscCall(PCBDDCSetupFETIDPMatContext(fetidpmat_ctx)); 2490 PetscCall(MatCreateShell(comm, fetidpmat_ctx->n, fetidpmat_ctx->n, fetidpmat_ctx->N, fetidpmat_ctx->N, fetidpmat_ctx, &newmat)); 2491 PetscCall(PetscObjectSetName((PetscObject)newmat, !fetidpmat_ctx->l2g_lambda_only ? "F" : "G")); 2492 PetscCall(MatShellSetOperation(newmat, MATOP_MULT, (void (*)(void))FETIDPMatMult)); 2493 PetscCall(MatShellSetOperation(newmat, MATOP_MULT_TRANSPOSE, (void (*)(void))FETIDPMatMultTranspose)); 2494 PetscCall(MatShellSetOperation(newmat, MATOP_DESTROY, (void (*)(void))PCBDDCDestroyFETIDPMat)); 2495 /* propagate MatOptions */ 2496 { 2497 PC_BDDC *pcbddc = (PC_BDDC *)fetidpmat_ctx->pc->data; 2498 PetscBool isset, issym; 2499 2500 PetscCall(MatIsSymmetricKnown(pc->mat, &isset, &issym)); 2501 if ((isset && issym) || pcbddc->symmetric_primal) PetscCall(MatSetOption(newmat, MAT_SYMMETRIC, PETSC_TRUE)); 2502 } 2503 PetscCall(MatSetOptionsPrefix(newmat, prefix)); 2504 PetscCall(MatAppendOptionsPrefix(newmat, "fetidp_")); 2505 PetscCall(MatSetUp(newmat)); 2506 /* FETI-DP preconditioner */ 2507 PetscCall(PCBDDCCreateFETIDPPCContext(pc, &fetidppc_ctx)); 2508 PetscCall(PCBDDCSetupFETIDPPCContext(newmat, fetidppc_ctx)); 2509 PetscCall(PCCreate(comm, &newpc)); 2510 PetscCall(PCSetOperators(newpc, newmat, newmat)); 2511 PetscCall(PCSetOptionsPrefix(newpc, prefix)); 2512 PetscCall(PCAppendOptionsPrefix(newpc, "fetidp_")); 2513 PetscCall(PCSetErrorIfFailure(newpc, pc->erroriffailure)); 2514 if (!fetidpmat_ctx->l2g_lambda_only) { /* standard FETI-DP */ 2515 PetscCall(PCSetType(newpc, PCSHELL)); 2516 PetscCall(PCShellSetName(newpc, "FETI-DP multipliers")); 2517 PetscCall(PCShellSetContext(newpc, fetidppc_ctx)); 2518 PetscCall(PCShellSetApply(newpc, FETIDPPCApply)); 2519 PetscCall(PCShellSetApplyTranspose(newpc, FETIDPPCApplyTranspose)); 2520 PetscCall(PCShellSetView(newpc, FETIDPPCView)); 2521 PetscCall(PCShellSetDestroy(newpc, PCBDDCDestroyFETIDPPC)); 2522 } else { /* saddle-point FETI-DP */ 2523 Mat M; 2524 PetscInt psize; 2525 PetscBool fake = PETSC_FALSE, isfieldsplit; 2526 2527 PetscCall(ISViewFromOptions(fetidpmat_ctx->lagrange, NULL, "-lag_view")); 2528 PetscCall(ISViewFromOptions(fetidpmat_ctx->pressure, NULL, "-press_view")); 2529 PetscCall(PetscObjectQuery((PetscObject)pc, "__KSPFETIDP_PPmat", (PetscObject *)&M)); 2530 PetscCall(PCSetType(newpc, PCFIELDSPLIT)); 2531 PetscCall(PCFieldSplitSetIS(newpc, "lag", fetidpmat_ctx->lagrange)); 2532 PetscCall(PCFieldSplitSetIS(newpc, "p", fetidpmat_ctx->pressure)); 2533 PetscCall(PCFieldSplitSetType(newpc, PC_COMPOSITE_SCHUR)); 2534 PetscCall(PCFieldSplitSetSchurFactType(newpc, PC_FIELDSPLIT_SCHUR_FACT_DIAG)); 2535 PetscCall(ISGetSize(fetidpmat_ctx->pressure, &psize)); 2536 if (psize != M->rmap->N) { 2537 Mat M2; 2538 PetscInt lpsize; 2539 2540 fake = PETSC_TRUE; 2541 PetscCall(ISGetLocalSize(fetidpmat_ctx->pressure, &lpsize)); 2542 PetscCall(MatCreate(comm, &M2)); 2543 PetscCall(MatSetType(M2, MATAIJ)); 2544 PetscCall(MatSetSizes(M2, lpsize, lpsize, psize, psize)); 2545 PetscCall(MatSetUp(M2)); 2546 PetscCall(MatAssemblyBegin(M2, MAT_FINAL_ASSEMBLY)); 2547 PetscCall(MatAssemblyEnd(M2, MAT_FINAL_ASSEMBLY)); 2548 PetscCall(PCFieldSplitSetSchurPre(newpc, PC_FIELDSPLIT_SCHUR_PRE_USER, M2)); 2549 PetscCall(MatDestroy(&M2)); 2550 } else { 2551 PetscCall(PCFieldSplitSetSchurPre(newpc, PC_FIELDSPLIT_SCHUR_PRE_USER, M)); 2552 } 2553 PetscCall(PCFieldSplitSetSchurScale(newpc, 1.0)); 2554 2555 /* we need to setfromoptions and setup here to access the blocks */ 2556 PetscCall(PCSetFromOptions(newpc)); 2557 PetscCall(PCSetUp(newpc)); 2558 2559 /* user may have changed the type (e.g. -fetidp_pc_type none) */ 2560 PetscCall(PetscObjectTypeCompare((PetscObject)newpc, PCFIELDSPLIT, &isfieldsplit)); 2561 if (isfieldsplit) { 2562 KSP *ksps; 2563 PC ppc, lagpc; 2564 PetscInt nn; 2565 PetscBool ismatis, matisok = PETSC_FALSE, check = PETSC_FALSE; 2566 2567 /* set the solver for the (0,0) block */ 2568 PetscCall(PCFieldSplitSchurGetSubKSP(newpc, &nn, &ksps)); 2569 if (!nn) { /* not of type PC_COMPOSITE_SCHUR */ 2570 PetscCall(PCFieldSplitGetSubKSP(newpc, &nn, &ksps)); 2571 if (!fake) { /* pass pmat to the pressure solver */ 2572 Mat F; 2573 2574 PetscCall(KSPGetOperators(ksps[1], &F, NULL)); 2575 PetscCall(KSPSetOperators(ksps[1], F, M)); 2576 } 2577 } else { 2578 PetscBool issym, isset; 2579 Mat S; 2580 2581 PetscCall(PCFieldSplitSchurGetS(newpc, &S)); 2582 PetscCall(MatIsSymmetricKnown(newmat, &isset, &issym)); 2583 if (isset) PetscCall(MatSetOption(S, MAT_SYMMETRIC, issym)); 2584 } 2585 PetscCall(KSPGetPC(ksps[0], &lagpc)); 2586 PetscCall(PCSetType(lagpc, PCSHELL)); 2587 PetscCall(PCShellSetName(lagpc, "FETI-DP multipliers")); 2588 PetscCall(PCShellSetContext(lagpc, fetidppc_ctx)); 2589 PetscCall(PCShellSetApply(lagpc, FETIDPPCApply)); 2590 PetscCall(PCShellSetApplyTranspose(lagpc, FETIDPPCApplyTranspose)); 2591 PetscCall(PCShellSetView(lagpc, FETIDPPCView)); 2592 PetscCall(PCShellSetDestroy(lagpc, PCBDDCDestroyFETIDPPC)); 2593 2594 /* Olof's idea: interface Schur complement preconditioner for the mass matrix */ 2595 PetscCall(KSPGetPC(ksps[1], &ppc)); 2596 if (fake) { 2597 BDDCIPC_ctx bddcipc_ctx; 2598 PetscContainer c; 2599 2600 matisok = PETSC_TRUE; 2601 2602 /* create inner BDDC solver */ 2603 PetscCall(PetscNew(&bddcipc_ctx)); 2604 PetscCall(PCCreate(comm, &bddcipc_ctx->bddc)); 2605 PetscCall(PCSetType(bddcipc_ctx->bddc, PCBDDC)); 2606 PetscCall(PCSetOperators(bddcipc_ctx->bddc, M, M)); 2607 PetscCall(PetscObjectQuery((PetscObject)pc, "__KSPFETIDP_pCSR", (PetscObject *)&c)); 2608 PetscCall(PetscObjectTypeCompare((PetscObject)M, MATIS, &ismatis)); 2609 if (c && ismatis) { 2610 Mat lM; 2611 PetscInt *csr, n; 2612 2613 PetscCall(MatISGetLocalMat(M, &lM)); 2614 PetscCall(MatGetSize(lM, &n, NULL)); 2615 PetscCall(PetscContainerGetPointer(c, (void **)&csr)); 2616 PetscCall(PCBDDCSetLocalAdjacencyGraph(bddcipc_ctx->bddc, n, csr, csr + (n + 1), PETSC_COPY_VALUES)); 2617 PetscCall(MatISRestoreLocalMat(M, &lM)); 2618 } 2619 PetscCall(PCSetOptionsPrefix(bddcipc_ctx->bddc, ((PetscObject)ksps[1])->prefix)); 2620 PetscCall(PCSetErrorIfFailure(bddcipc_ctx->bddc, pc->erroriffailure)); 2621 PetscCall(PCSetFromOptions(bddcipc_ctx->bddc)); 2622 2623 /* wrap the interface application */ 2624 PetscCall(PCSetType(ppc, PCSHELL)); 2625 PetscCall(PCShellSetName(ppc, "FETI-DP pressure")); 2626 PetscCall(PCShellSetContext(ppc, bddcipc_ctx)); 2627 PetscCall(PCShellSetSetUp(ppc, PCSetUp_BDDCIPC)); 2628 PetscCall(PCShellSetApply(ppc, PCApply_BDDCIPC)); 2629 PetscCall(PCShellSetApplyTranspose(ppc, PCApplyTranspose_BDDCIPC)); 2630 PetscCall(PCShellSetView(ppc, PCView_BDDCIPC)); 2631 PetscCall(PCShellSetDestroy(ppc, PCDestroy_BDDCIPC)); 2632 } 2633 2634 /* determine if we need to assemble M to construct a preconditioner */ 2635 if (!matisok) { 2636 PetscCall(PetscObjectTypeCompare((PetscObject)M, MATIS, &ismatis)); 2637 PetscCall(PetscObjectTypeCompareAny((PetscObject)ppc, &matisok, PCBDDC, PCJACOBI, PCNONE, PCMG, "")); 2638 if (ismatis && !matisok) PetscCall(MatConvert(M, MATAIJ, MAT_INPLACE_MATRIX, &M)); 2639 } 2640 2641 /* run the subproblems to check convergence */ 2642 PetscCall(PetscOptionsGetBool(NULL, ((PetscObject)newmat)->prefix, "-check_saddlepoint", &check, NULL)); 2643 if (check) { 2644 PetscInt i; 2645 2646 for (i = 0; i < nn; i++) { 2647 KSP kspC; 2648 PC pc; 2649 Mat F, pF; 2650 Vec x, y; 2651 PetscBool isschur, prec = PETSC_TRUE; 2652 2653 PetscCall(KSPCreate(PetscObjectComm((PetscObject)ksps[i]), &kspC)); 2654 PetscCall(KSPSetOptionsPrefix(kspC, ((PetscObject)ksps[i])->prefix)); 2655 PetscCall(KSPAppendOptionsPrefix(kspC, "check_")); 2656 PetscCall(KSPGetOperators(ksps[i], &F, &pF)); 2657 PetscCall(PetscObjectTypeCompare((PetscObject)F, MATSCHURCOMPLEMENT, &isschur)); 2658 if (isschur) { 2659 KSP kspS, kspS2; 2660 Mat A00, pA00, A10, A01, A11; 2661 char prefix[256]; 2662 2663 PetscCall(MatSchurComplementGetKSP(F, &kspS)); 2664 PetscCall(MatSchurComplementGetSubMatrices(F, &A00, &pA00, &A01, &A10, &A11)); 2665 PetscCall(MatCreateSchurComplement(A00, pA00, A01, A10, A11, &F)); 2666 PetscCall(MatSchurComplementGetKSP(F, &kspS2)); 2667 PetscCall(PetscSNPrintf(prefix, sizeof(prefix), "%sschur_", ((PetscObject)kspC)->prefix)); 2668 PetscCall(KSPSetOptionsPrefix(kspS2, prefix)); 2669 PetscCall(KSPGetPC(kspS2, &pc)); 2670 PetscCall(PCSetType(pc, PCKSP)); 2671 PetscCall(PCKSPSetKSP(pc, kspS)); 2672 PetscCall(KSPSetFromOptions(kspS2)); 2673 PetscCall(KSPGetPC(kspS2, &pc)); 2674 PetscCall(PCSetUseAmat(pc, PETSC_TRUE)); 2675 } else { 2676 PetscCall(PetscObjectReference((PetscObject)F)); 2677 } 2678 PetscCall(KSPSetFromOptions(kspC)); 2679 PetscCall(PetscOptionsGetBool(NULL, ((PetscObject)kspC)->prefix, "-preconditioned", &prec, NULL)); 2680 if (prec) { 2681 PetscCall(KSPGetPC(ksps[i], &pc)); 2682 PetscCall(KSPSetPC(kspC, pc)); 2683 } 2684 PetscCall(KSPSetOperators(kspC, F, pF)); 2685 PetscCall(MatCreateVecs(F, &x, &y)); 2686 PetscCall(VecSetRandom(x, NULL)); 2687 PetscCall(MatMult(F, x, y)); 2688 PetscCall(KSPSolve(kspC, y, x)); 2689 PetscCall(KSPCheckSolve(kspC, pc, x)); 2690 PetscCall(KSPDestroy(&kspC)); 2691 PetscCall(MatDestroy(&F)); 2692 PetscCall(VecDestroy(&x)); 2693 PetscCall(VecDestroy(&y)); 2694 } 2695 } 2696 PetscCall(PetscFree(ksps)); 2697 } 2698 } 2699 /* return pointers for objects created */ 2700 *fetidp_mat = newmat; 2701 *fetidp_pc = newpc; 2702 PetscFunctionReturn(0); 2703 } 2704 2705 /*@C 2706 PCBDDCCreateFETIDPOperators - Create FETI-DP operators 2707 2708 Collective 2709 2710 Input Parameters: 2711 + pc - the BDDC preconditioning context (setup should have been called before) 2712 . fully_redundant - true for a fully redundant set of Lagrange multipliers 2713 - prefix - optional options database prefix for the objects to be created (can be NULL) 2714 2715 Output Parameters: 2716 + fetidp_mat - shell FETI-DP matrix object 2717 - fetidp_pc - shell Dirichlet preconditioner for FETI-DP matrix 2718 2719 Level: developer 2720 2721 Note: 2722 Currently the only operations provided for FETI-DP matrix are `MatMult()` and `MatMultTranspose()` 2723 2724 .seealso: `PCBDDC`, `PCBDDCMatFETIDPGetRHS()`, `PCBDDCMatFETIDPGetSolution()` 2725 @*/ 2726 PetscErrorCode PCBDDCCreateFETIDPOperators(PC pc, PetscBool fully_redundant, const char *prefix, Mat *fetidp_mat, PC *fetidp_pc) 2727 { 2728 PetscFunctionBegin; 2729 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 2730 if (pc->setupcalled) { 2731 PetscUseMethod(pc, "PCBDDCCreateFETIDPOperators_C", (PC, PetscBool, const char *, Mat *, PC *), (pc, fully_redundant, prefix, fetidp_mat, fetidp_pc)); 2732 } else SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "You must call PCSetup_BDDC() first"); 2733 PetscFunctionReturn(0); 2734 } 2735 2736 /*MC 2737 PCBDDC - Balancing Domain Decomposition by Constraints preconditioners 2738 2739 Requires `MATIS` matrices (Pmat) with local matrices (inside the `MATIS`) of type `MATSEQAIJ`, `MATSEQBAIJ` or `MATSEQSBAIJ` 2740 2741 It also works with unsymmetric and indefinite problems. 2742 2743 Unlike 'conventional' interface preconditioners, `PCBDDC` iterates over all degrees of freedom, not just those on the interface. This allows the use 2744 of approximate solvers on the subdomains. 2745 2746 Approximate local solvers are automatically adapted (see [1]) if the user has attached a nullspace object to the subdomain matrices, and informed 2747 `PCBDDC` of using approximate solvers (via the command line). 2748 2749 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. 2750 The latter can be customized by using `PCBDDCSetLocalAdjacencyGraph()` 2751 2752 Additional information on dofs can be provided by using `PCBDDCSetDofsSplitting()`, `PCBDDCSetDirichletBoundaries()`, `PCBDDCSetNeumannBoundaries()`, and 2753 `PCBDDCSetPrimalVerticesIS()` and their local counterparts. 2754 2755 Constraints can be customized by attaching a `MatNullSpace` object to the `MATIS` matrix via `MatSetNearNullSpace()`. Non-singular modes are retained via SVD. 2756 2757 Change of basis is performed similarly to [2] when requested. When more than one constraint is present on a single connected component 2758 (i.e. an edge or a face), a robust method based on local QR factorizations is used. 2759 User defined change of basis can be passed to `PCBDDC` with `PCBDDCSetChangeOfBasisMat()` 2760 2761 The PETSc implementation also supports multilevel `PCBDDC` [3]. Coarse grids are partitioned using a `MatPartitioning` object. 2762 2763 Adaptive selection of primal constraints [4] is supported for SPD systems with high-contrast in the coefficients if MUMPS or MKL_PARDISO are present. 2764 Future versions of the code will also consider using PASTIX. 2765 2766 An experimental interface to the FETI-DP method is available. FETI-DP operators could be created using `PCBDDCCreateFETIDPOperators()`. 2767 A stand-alone class for the FETI-DP method will be provided in the next releases. 2768 2769 Options Database Keys: 2770 + -pc_bddc_use_vertices <true> - use or not vertices in primal space 2771 . -pc_bddc_use_edges <true> - use or not edges in primal space 2772 . -pc_bddc_use_faces <false> - use or not faces in primal space 2773 . -pc_bddc_symmetric <true> - symmetric computation of primal basis functions. Specify false for unsymmetric problems 2774 . -pc_bddc_use_change_of_basis <false> - use change of basis approach (on edges only) 2775 . -pc_bddc_use_change_on_faces <false> - use change of basis approach on faces if change of basis has been requested 2776 . -pc_bddc_switch_static <false> - switches from M_2 (default) to M_3 operator (see reference article [1]) 2777 . -pc_bddc_levels <0> - maximum number of levels for multilevel 2778 . -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) 2779 . -pc_bddc_coarse_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) 2780 . -pc_bddc_use_deluxe_scaling <false> - use deluxe scaling 2781 . -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) 2782 . -pc_bddc_adaptive_threshold <0.0> - when a value different than zero is specified, adaptive selection of constraints is performed on edges and faces (requires deluxe scaling and MUMPS or MKL_PARDISO installed) 2783 - -pc_bddc_check_level <0> - set verbosity level of debugging output 2784 2785 Options for Dirichlet, Neumann or coarse solver can be set using the appropriate options prefix 2786 .vb 2787 -pc_bddc_dirichlet_ 2788 -pc_bddc_neumann_ 2789 -pc_bddc_coarse_ 2790 .ve 2791 e.g. -pc_bddc_dirichlet_ksp_type richardson -pc_bddc_dirichlet_pc_type gamg. `PCBDDC` uses by default `KSPPREONLY` and `PCLU`. 2792 2793 When using a multilevel approach, solvers' options at the N-th level (N > 1) can be specified using the options prefix 2794 .vb 2795 -pc_bddc_dirichlet_lN_ 2796 -pc_bddc_neumann_lN_ 2797 -pc_bddc_coarse_lN_ 2798 .ve 2799 Note that level number ranges from the finest (0) to the coarsest (N). 2800 In order to specify options for the `PCBDDC` operators at the coarser levels (and not for the solvers), prepend -pc_bddc_coarse_ or -pc_bddc_coarse_l 2801 to the option, e.g. 2802 .vb 2803 -pc_bddc_coarse_pc_bddc_adaptive_threshold 5 -pc_bddc_coarse_l1_pc_bddc_redistribute 3 2804 .ve 2805 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 2806 2807 References: 2808 + * - C. R. Dohrmann. "An approximate BDDC preconditioner", Numerical Linear Algebra with Applications Volume 14, Issue 2, pages 149--168, March 2007 2809 . * - A. Klawonn and O. B. Widlund. "Dual-Primal FETI Methods for Linear Elasticity", Communications on Pure and Applied Mathematics Volume 59, Issue 11, pages 1523--1572, November 2006 2810 . * - J. Mandel, B. Sousedik, C. R. Dohrmann. "Multispace and Multilevel BDDC", Computing Volume 83, Issue 2--3, pages 55--85, November 2008 2811 - * - 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 2812 2813 Level: intermediate 2814 2815 Contributed by Stefano Zampini 2816 2817 .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `MATIS`, `PCLU`, `PGGAMG`, `PC`, `PCBDDCSetLocalAdjacencyGraph()`, `PCBDDCSetDofsSplitting()`, 2818 `PCBDDCSetDirichletBoundaries()`, `PCBDDCSetNeumannBoundaries()`, `PCBDDCSetPrimalVerticesIS()`, `MatNullSpace`, `MatSetNearNullSpace()`, 2819 `PCBDDCSetChangeOfBasisMat()`, `PCBDDCCreateFETIDPOperators()`, `PCNN` 2820 M*/ 2821 2822 PETSC_EXTERN PetscErrorCode PCCreate_BDDC(PC pc) 2823 { 2824 PC_BDDC *pcbddc; 2825 2826 PetscFunctionBegin; 2827 PetscCall(PetscNew(&pcbddc)); 2828 pc->data = pcbddc; 2829 2830 /* create PCIS data structure */ 2831 PetscCall(PCISCreate(pc)); 2832 2833 /* create local graph structure */ 2834 PetscCall(PCBDDCGraphCreate(&pcbddc->mat_graph)); 2835 2836 /* BDDC nonzero defaults */ 2837 pcbddc->use_nnsp = PETSC_TRUE; 2838 pcbddc->use_local_adj = PETSC_TRUE; 2839 pcbddc->use_vertices = PETSC_TRUE; 2840 pcbddc->use_edges = PETSC_TRUE; 2841 pcbddc->symmetric_primal = PETSC_TRUE; 2842 pcbddc->vertex_size = 1; 2843 pcbddc->recompute_topography = PETSC_TRUE; 2844 pcbddc->coarse_size = -1; 2845 pcbddc->use_exact_dirichlet_trick = PETSC_TRUE; 2846 pcbddc->coarsening_ratio = 8; 2847 pcbddc->coarse_eqs_per_proc = 1; 2848 pcbddc->benign_compute_correction = PETSC_TRUE; 2849 pcbddc->nedfield = -1; 2850 pcbddc->nedglobal = PETSC_TRUE; 2851 pcbddc->graphmaxcount = PETSC_MAX_INT; 2852 pcbddc->sub_schurs_layers = -1; 2853 pcbddc->adaptive_threshold[0] = 0.0; 2854 pcbddc->adaptive_threshold[1] = 0.0; 2855 2856 /* function pointers */ 2857 pc->ops->apply = PCApply_BDDC; 2858 pc->ops->applytranspose = PCApplyTranspose_BDDC; 2859 pc->ops->setup = PCSetUp_BDDC; 2860 pc->ops->destroy = PCDestroy_BDDC; 2861 pc->ops->setfromoptions = PCSetFromOptions_BDDC; 2862 pc->ops->view = PCView_BDDC; 2863 pc->ops->applyrichardson = NULL; 2864 pc->ops->applysymmetricleft = NULL; 2865 pc->ops->applysymmetricright = NULL; 2866 pc->ops->presolve = PCPreSolve_BDDC; 2867 pc->ops->postsolve = PCPostSolve_BDDC; 2868 pc->ops->reset = PCReset_BDDC; 2869 2870 /* composing function */ 2871 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBDDCSetDiscreteGradient_C", PCBDDCSetDiscreteGradient_BDDC)); 2872 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBDDCSetDivergenceMat_C", PCBDDCSetDivergenceMat_BDDC)); 2873 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBDDCSetChangeOfBasisMat_C", PCBDDCSetChangeOfBasisMat_BDDC)); 2874 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBDDCSetPrimalVerticesLocalIS_C", PCBDDCSetPrimalVerticesLocalIS_BDDC)); 2875 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBDDCSetPrimalVerticesIS_C", PCBDDCSetPrimalVerticesIS_BDDC)); 2876 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBDDCGetPrimalVerticesLocalIS_C", PCBDDCGetPrimalVerticesLocalIS_BDDC)); 2877 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBDDCGetPrimalVerticesIS_C", PCBDDCGetPrimalVerticesIS_BDDC)); 2878 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBDDCSetCoarseningRatio_C", PCBDDCSetCoarseningRatio_BDDC)); 2879 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBDDCSetLevel_C", PCBDDCSetLevel_BDDC)); 2880 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBDDCSetUseExactDirichlet_C", PCBDDCSetUseExactDirichlet_BDDC)); 2881 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBDDCSetLevels_C", PCBDDCSetLevels_BDDC)); 2882 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBDDCSetDirichletBoundaries_C", PCBDDCSetDirichletBoundaries_BDDC)); 2883 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBDDCSetDirichletBoundariesLocal_C", PCBDDCSetDirichletBoundariesLocal_BDDC)); 2884 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBDDCSetNeumannBoundaries_C", PCBDDCSetNeumannBoundaries_BDDC)); 2885 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBDDCSetNeumannBoundariesLocal_C", PCBDDCSetNeumannBoundariesLocal_BDDC)); 2886 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBDDCGetDirichletBoundaries_C", PCBDDCGetDirichletBoundaries_BDDC)); 2887 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBDDCGetDirichletBoundariesLocal_C", PCBDDCGetDirichletBoundariesLocal_BDDC)); 2888 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBDDCGetNeumannBoundaries_C", PCBDDCGetNeumannBoundaries_BDDC)); 2889 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBDDCGetNeumannBoundariesLocal_C", PCBDDCGetNeumannBoundariesLocal_BDDC)); 2890 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBDDCSetDofsSplitting_C", PCBDDCSetDofsSplitting_BDDC)); 2891 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBDDCSetDofsSplittingLocal_C", PCBDDCSetDofsSplittingLocal_BDDC)); 2892 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBDDCSetLocalAdjacencyGraph_C", PCBDDCSetLocalAdjacencyGraph_BDDC)); 2893 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBDDCCreateFETIDPOperators_C", PCBDDCCreateFETIDPOperators_BDDC)); 2894 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBDDCMatFETIDPGetRHS_C", PCBDDCMatFETIDPGetRHS_BDDC)); 2895 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBDDCMatFETIDPGetSolution_C", PCBDDCMatFETIDPGetSolution_BDDC)); 2896 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCPreSolveChangeRHS_C", PCPreSolveChangeRHS_BDDC)); 2897 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSetCoordinates_C", PCSetCoordinates_BDDC)); 2898 PetscFunctionReturn(0); 2899 } 2900 2901 /*@C 2902 PCBDDCInitializePackage - This function initializes everything in the `PCBDDC` package. It is called 2903 from `PCInitializePackage()`. 2904 2905 Level: developer 2906 2907 .seealso: `PetscInitialize()`, `PCBDDCFinalizePackage()` 2908 @*/ 2909 PetscErrorCode PCBDDCInitializePackage(void) 2910 { 2911 int i; 2912 2913 PetscFunctionBegin; 2914 if (PCBDDCPackageInitialized) PetscFunctionReturn(0); 2915 PCBDDCPackageInitialized = PETSC_TRUE; 2916 PetscCall(PetscRegisterFinalize(PCBDDCFinalizePackage)); 2917 2918 /* general events */ 2919 PetscCall(PetscLogEventRegister("PCBDDCTopo", PC_CLASSID, &PC_BDDC_Topology[0])); 2920 PetscCall(PetscLogEventRegister("PCBDDCLKSP", PC_CLASSID, &PC_BDDC_LocalSolvers[0])); 2921 PetscCall(PetscLogEventRegister("PCBDDCLWor", PC_CLASSID, &PC_BDDC_LocalWork[0])); 2922 PetscCall(PetscLogEventRegister("PCBDDCCorr", PC_CLASSID, &PC_BDDC_CorrectionSetUp[0])); 2923 PetscCall(PetscLogEventRegister("PCBDDCASet", PC_CLASSID, &PC_BDDC_ApproxSetUp[0])); 2924 PetscCall(PetscLogEventRegister("PCBDDCAApp", PC_CLASSID, &PC_BDDC_ApproxApply[0])); 2925 PetscCall(PetscLogEventRegister("PCBDDCCSet", PC_CLASSID, &PC_BDDC_CoarseSetUp[0])); 2926 PetscCall(PetscLogEventRegister("PCBDDCCKSP", PC_CLASSID, &PC_BDDC_CoarseSolver[0])); 2927 PetscCall(PetscLogEventRegister("PCBDDCAdap", PC_CLASSID, &PC_BDDC_AdaptiveSetUp[0])); 2928 PetscCall(PetscLogEventRegister("PCBDDCScal", PC_CLASSID, &PC_BDDC_Scaling[0])); 2929 PetscCall(PetscLogEventRegister("PCBDDCSchr", PC_CLASSID, &PC_BDDC_Schurs[0])); 2930 PetscCall(PetscLogEventRegister("PCBDDCDirS", PC_CLASSID, &PC_BDDC_Solves[0][0])); 2931 PetscCall(PetscLogEventRegister("PCBDDCNeuS", PC_CLASSID, &PC_BDDC_Solves[0][1])); 2932 PetscCall(PetscLogEventRegister("PCBDDCCoaS", PC_CLASSID, &PC_BDDC_Solves[0][2])); 2933 for (i = 1; i < PETSC_PCBDDC_MAXLEVELS; i++) { 2934 char ename[32]; 2935 2936 PetscCall(PetscSNPrintf(ename, sizeof(ename), "PCBDDCTopo l%02d", i)); 2937 PetscCall(PetscLogEventRegister(ename, PC_CLASSID, &PC_BDDC_Topology[i])); 2938 PetscCall(PetscSNPrintf(ename, sizeof(ename), "PCBDDCLKSP l%02d", i)); 2939 PetscCall(PetscLogEventRegister(ename, PC_CLASSID, &PC_BDDC_LocalSolvers[i])); 2940 PetscCall(PetscSNPrintf(ename, sizeof(ename), "PCBDDCLWor l%02d", i)); 2941 PetscCall(PetscLogEventRegister(ename, PC_CLASSID, &PC_BDDC_LocalWork[i])); 2942 PetscCall(PetscSNPrintf(ename, sizeof(ename), "PCBDDCCorr l%02d", i)); 2943 PetscCall(PetscLogEventRegister(ename, PC_CLASSID, &PC_BDDC_CorrectionSetUp[i])); 2944 PetscCall(PetscSNPrintf(ename, sizeof(ename), "PCBDDCASet l%02d", i)); 2945 PetscCall(PetscLogEventRegister(ename, PC_CLASSID, &PC_BDDC_ApproxSetUp[i])); 2946 PetscCall(PetscSNPrintf(ename, sizeof(ename), "PCBDDCAApp l%02d", i)); 2947 PetscCall(PetscLogEventRegister(ename, PC_CLASSID, &PC_BDDC_ApproxApply[i])); 2948 PetscCall(PetscSNPrintf(ename, sizeof(ename), "PCBDDCCSet l%02d", i)); 2949 PetscCall(PetscLogEventRegister(ename, PC_CLASSID, &PC_BDDC_CoarseSetUp[i])); 2950 PetscCall(PetscSNPrintf(ename, sizeof(ename), "PCBDDCCKSP l%02d", i)); 2951 PetscCall(PetscLogEventRegister(ename, PC_CLASSID, &PC_BDDC_CoarseSolver[i])); 2952 PetscCall(PetscSNPrintf(ename, sizeof(ename), "PCBDDCAdap l%02d", i)); 2953 PetscCall(PetscLogEventRegister(ename, PC_CLASSID, &PC_BDDC_AdaptiveSetUp[i])); 2954 PetscCall(PetscSNPrintf(ename, sizeof(ename), "PCBDDCScal l%02d", i)); 2955 PetscCall(PetscLogEventRegister(ename, PC_CLASSID, &PC_BDDC_Scaling[i])); 2956 PetscCall(PetscSNPrintf(ename, sizeof(ename), "PCBDDCSchr l%02d", i)); 2957 PetscCall(PetscLogEventRegister(ename, PC_CLASSID, &PC_BDDC_Schurs[i])); 2958 PetscCall(PetscSNPrintf(ename, sizeof(ename), "PCBDDCDirS l%02d", i)); 2959 PetscCall(PetscLogEventRegister(ename, PC_CLASSID, &PC_BDDC_Solves[i][0])); 2960 PetscCall(PetscSNPrintf(ename, sizeof(ename), "PCBDDCNeuS l%02d", i)); 2961 PetscCall(PetscLogEventRegister(ename, PC_CLASSID, &PC_BDDC_Solves[i][1])); 2962 PetscCall(PetscSNPrintf(ename, sizeof(ename), "PCBDDCCoaS l%02d", i)); 2963 PetscCall(PetscLogEventRegister(ename, PC_CLASSID, &PC_BDDC_Solves[i][2])); 2964 } 2965 PetscFunctionReturn(0); 2966 } 2967 2968 /*@C 2969 PCBDDCFinalizePackage - This function frees everything from the `PCBDDC` package. It is 2970 called from `PetscFinalize()` automatically. 2971 2972 Level: developer 2973 2974 .seealso: `PetscFinalize()`, `PCBDDCInitializePackage()` 2975 @*/ 2976 PetscErrorCode PCBDDCFinalizePackage(void) 2977 { 2978 PetscFunctionBegin; 2979 PCBDDCPackageInitialized = PETSC_FALSE; 2980 PetscFunctionReturn(0); 2981 } 2982