xref: /petsc/include/petsc/private/pcbddcstructsimpl.h (revision f0d2bb25642116b571c72fa9797d32d784b1f106)
1 #pragma once
2 
3 #include <petscksp.h>
4 #include <petscbt.h>
5 
6 /* special marks for interface graph: they cannot be enums
7    since PCBDDCGRAPH_SPECIAL_MARK ranges from -4 to -max_int */
8 #define PCBDDCGRAPH_NEUMANN_MARK        -1
9 #define PCBDDCGRAPH_DIRICHLET_MARK      -2
10 #define PCBDDCGRAPH_LOCAL_PERIODIC_MARK -3
11 #define PCBDDCGRAPH_SPECIAL_MARK        -4
12 
13 /* Metadata information on node */
14 typedef struct {
15   PetscBool touched;
16   PetscInt  subset;
17   PetscInt  which_dof;
18   PetscInt  special_dof;
19   PetscInt  local_sub;
20   PetscInt  count;
21   PetscInt *neighbours_set;
22   PetscBool shared;
23   PetscInt  local_groups_count;
24   PetscInt *local_groups;
25 } PCBDDCGraphNode;
26 
27 /* Data structure for local graph partitioning */
28 struct _PCBDDCGraph {
29   PetscBool setupcalled;
30   /* graph information */
31   ISLocalToGlobalMapping l2gmap;
32   PetscInt               nvtxs;
33   PetscInt               nvtxs_global;
34   PCBDDCGraphNode       *nodes;
35   PetscInt               custom_minimal_size;
36   PetscBool              twodim;
37   PetscBool              twodimset;
38   PetscBool              has_dirichlet;
39   PetscBool              multi_element;
40   IS                     dirdofs;
41   IS                     dirdofsB;
42   PetscBool              seq_graph;
43   PetscInt               maxcount;
44   /* data for connected components */
45   PetscInt  ncc;
46   PetscInt *cptr;
47   PetscInt *queue;
48   PetscBool queue_sorted;
49   /* data for interface subsets */
50   PetscInt   n_subsets;
51   PetscInt  *subset_size;
52   PetscInt **subset_idxs;
53   PetscInt  *subset_ncc;
54   PetscInt  *subset_ref_node;
55   PetscInt  *gsubset_size;
56   PetscInt  *interface_ref_rsize;
57   PetscSF    interface_ref_sf;
58   PetscSF    interface_subset_sf;
59   /* placeholders for connectivity relation between dofs */
60   PetscInt  nvtxs_csr;
61   PetscInt *xadj;
62   PetscInt *adjncy;
63   PetscBool freecsr;
64   /* data for local subdomains (if any have been detected)
65      these are not intended to be exposed */
66   PetscInt  n_local_subs;
67   PetscInt *local_subs;
68   /* coordinates (for corner detection) */
69   PetscBool  active_coords;
70   PetscBool  cloc;
71   PetscInt   cdim, cnloc;
72   PetscReal *coords;
73 };
74 typedef struct _PCBDDCGraph *PCBDDCGraph;
75 
76 struct _PCBDDCGraphCandidates {
77   PetscInt nfc, nec;
78   IS      *Faces, *Edges, Vertices;
79 };
80 typedef struct _PCBDDCGraphCandidates *PCBDDCGraphCandidates;
81 
82 /* Wrap to MatFactor solver in Schur complement mode. Provides
83    - standalone solver for interior variables
84    - forward and backward substitutions for correction solver
85 */
86 /* It assumes that interior variables are a contiguous set starting from 0 */
87 struct _PCBDDCReuseSolvers {
88   /* the factored matrix obtained from MatGetFactor(...,solver_package,...) */
89   Mat F;
90   /* placeholders for the solution and rhs on the whole set of dofs of A (size local_dofs - local_vertices)*/
91   Vec sol;
92   Vec rhs;
93   /* */
94   PetscBool has_vertices;
95   /* shell PCs to handle interior/correction solvers */
96   PC interior_solver;
97   PC correction_solver;
98   IS is_R;
99   /* objects to handle Schur complement solution */
100   Vec        rhs_B;
101   Vec        sol_B;
102   IS         is_B;
103   VecScatter correction_scatter_B;
104   /* handle benign trick without change of basis on pressures */
105   PetscInt     benign_n;
106   IS          *benign_zerodiag_subs;
107   PetscScalar *benign_save_vals;
108   Mat          benign_csAIB;
109   Mat          benign_AIIm1ones;
110   Vec          benign_corr_work;
111   Vec          benign_dummy_schur_vec;
112 };
113 typedef struct _PCBDDCReuseSolvers *PCBDDCReuseSolvers;
114 
115 /* structure to handle Schur complements on subsets */
116 struct _PCBDDCSubSchurs {
117   /* local Neumann matrix */
118   Mat A;
119   /* local Schur complement */
120   Mat S;
121   /* index sets */
122   IS is_I;
123   IS is_B;
124   /* whether Schur complements are explicitly computed with or not */
125   PetscBool schur_explicit;
126   /* BDDC or GDSW */
127   PetscBool gdsw;
128   /* matrices contained explicit schur complements cat together */
129   /* note that AIJ format is used but the values are inserted as in column major ordering */
130   Mat S_Ej_all;
131   Mat sum_S_Ej_all;
132   Mat sum_S_Ej_inv_all;
133   Mat sum_S_Ej_tilda_all;
134   IS  is_Ej_all;
135   IS  is_vertices;
136   IS  is_dir;
137   /* l2g maps */
138   ISLocalToGlobalMapping l2gmap;
139   ISLocalToGlobalMapping BtoNmap;
140   /* number of local subproblems */
141   PetscInt n_subs;
142   /* connected components */
143   IS     *is_subs;
144   PetscBT is_edge;
145   /* mat factor */
146   char          mat_solver_type[64];
147   MatFactorType mat_factor_type;
148   PetscBool     is_symmetric;
149   PetscBool     is_hermitian;
150   PetscBool     is_posdef;
151   /* data structure to reuse MatFactor with Schur solver */
152   PCBDDCReuseSolvers reuse_solver;
153   /* change of variables */
154   KSP      *change;
155   IS       *change_primal_sub;
156   PetscBool change_with_qr;
157   /* prefix */
158   char *prefix;
159   /* */
160   PetscBool restrict_comm;
161   /* debug */
162   PetscBool debug;
163   /* hook to PCBDDCGraph */
164   PCBDDCGraph graph;
165 };
166 typedef struct _PCBDDCSubSchurs *PCBDDCSubSchurs;
167 
168 /* Structure for deluxe scaling */
169 struct _PCBDDCDeluxeScaling {
170   /* simple scaling on selected dofs (i.e. primal vertices and nodes on interface dirichlet boundaries) */
171   PetscInt  n_simple;
172   PetscInt *idx_simple_B;
173   /* handle deluxe problems  */
174   PetscInt     seq_n;
175   PetscScalar *workspace;
176   VecScatter  *seq_scctx;
177   Vec         *seq_work1;
178   Vec         *seq_work2;
179   Mat         *seq_mat;
180   Mat         *seq_mat_inv_sum;
181   KSP         *change;
182   PetscBool    change_with_qr;
183 };
184 typedef struct _PCBDDCDeluxeScaling *PCBDDCDeluxeScaling;
185 
186 /* inexact solvers with nullspace correction */
187 struct _NullSpaceCorrection_ctx {
188   Mat           basis_mat;
189   Mat           inv_smat;
190   PC            local_pc;
191   Vec          *fw;
192   Vec          *sw;
193   PetscScalar   scale;
194   PetscLogEvent evapply;
195   PetscBool     symm;
196 };
197 typedef struct _NullSpaceCorrection_ctx *NullSpaceCorrection_ctx;
198 
199 /* MatShell context for benign mat mults */
200 struct _PCBDDCBenignMatMult_ctx {
201   Mat          A;
202   PetscInt     benign_n;
203   IS          *benign_zerodiag_subs;
204   PetscScalar *work;
205   PetscBool    apply_left;
206   PetscBool    apply_right;
207   PetscBool    apply_p0;
208   PetscBool    free;
209 };
210 typedef struct _PCBDDCBenignMatMult_ctx *PCBDDCBenignMatMult_ctx;
211 
212 /* feti-dp mat */
213 struct _FETIDPMat_ctx {
214   PetscInt   n;        /* local number of rows */
215   PetscInt   N;        /* global number of rows */
216   PetscInt   n_lambda; /* global number of multipliers */
217   Vec        lambda_local;
218   Vec        temp_solution_B;
219   Vec        temp_solution_D;
220   Mat        B_delta;
221   Mat        B_Ddelta;
222   PetscBool  deluxe_nonred;
223   VecScatter l2g_lambda;
224   PC         pc;
225   PetscBool  fully_redundant;
226   /* saddle point */
227   VecScatter l2g_lambda_only;
228   Mat        B_BB;
229   Mat        B_BI;
230   Mat        Bt_BB;
231   Mat        Bt_BI;
232   Mat        C;
233   VecScatter l2g_p;
234   VecScatter g2g_p;
235   Vec        vP;
236   Vec        xPg;
237   Vec        yPg;
238   Vec        rhs_flip;
239   IS         lP_I;
240   IS         lP_B;
241   IS         pressure;
242   IS         lagrange;
243 };
244 typedef struct _FETIDPMat_ctx *FETIDPMat_ctx;
245 
246 /* feti-dp preconditioner */
247 struct _FETIDPPC_ctx {
248   Mat        S_j;
249   Vec        lambda_local;
250   Mat        B_Ddelta;
251   VecScatter l2g_lambda;
252   PC         pc;
253   /* saddle point */
254   Vec xPg;
255   Vec yPg;
256 };
257 typedef struct _FETIDPPC_ctx *FETIDPPC_ctx;
258 
259 struct _BDdelta_DN {
260   Mat BD;
261   KSP kBD;
262   Vec work;
263 };
264 typedef struct _BDdelta_DN *BDdelta_DN;
265 
266 /* Schur interface preconditioner */
267 struct _BDDCIPC_ctx {
268   VecScatter g2l;
269   PC         bddc;
270 };
271 typedef struct _BDDCIPC_ctx *BDDCIPC_ctx;
272