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