xref: /petsc/src/ksp/ksp/impls/fetidp/fetidp.c (revision 4977266fffaaedcd41e5a67923646569c4e329aa)
1 #include <petsc/private/kspimpl.h> /*I <petscksp.h> I*/
2 #include <petsc/private/pcbddcimpl.h>
3 #include <petsc/private/pcbddcprivateimpl.h>
4 #include <petscdm.h>
5 
6 static PetscBool  cited       = PETSC_FALSE;
7 static PetscBool  cited2      = 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"
20                                 "@article{ZampiniDualPrimal,\n"
21                                 "author = {Stefano Zampini},\n"
22                                 "title = {{D}ual-{P}rimal methods for the cardiac {B}idomain model},\n"
23                                 "volume = {24},\n"
24                                 "number = {04},\n"
25                                 "pages = {667-696},\n"
26                                 "year = {2014},\n"
27                                 "doi = {10.1142/S0218202513500632},\n"
28                                 "URL = {https://www.worldscientific.com/doi/abs/10.1142/S0218202513500632},\n"
29                                 "eprint = {https://www.worldscientific.com/doi/pdf/10.1142/S0218202513500632}\n"
30                                 "}\n";
31 static const char citation2[] = "@article{li2013nonoverlapping,\n"
32                                 "title={A nonoverlapping domain decomposition method for incompressible Stokes equations with continuous pressures},\n"
33                                 "author={Li, Jing and Tu, Xuemin},\n"
34                                 "journal={SIAM Journal on Numerical Analysis},\n"
35                                 "volume={51},\n"
36                                 "number={2},\n"
37                                 "pages={1235--1253},\n"
38                                 "year={2013},\n"
39                                 "publisher={Society for Industrial and Applied Mathematics}\n"
40                                 "}\n";
41 
42 /*
43     This file implements the FETI-DP method in PETSc as part of KSP.
44 */
45 typedef struct {
46   KSP parentksp;
47 } KSP_FETIDPMon;
48 
49 typedef struct {
50   KSP              innerksp;        /* the KSP for the Lagrange multipliers */
51   PC               innerbddc;       /* the inner BDDC object */
52   PetscBool        fully_redundant; /* true for using a fully redundant set of multipliers */
53   PetscBool        userbddc;        /* true if the user provided the PCBDDC object */
54   PetscBool        saddlepoint;     /* support for saddle point problems */
55   IS               pP;              /* index set for pressure variables */
56   Vec              rhs_flip;        /* see KSPFETIDPSetUpOperators */
57   KSP_FETIDPMon   *monctx;          /* monitor context, used to pass user defined monitors
58                                         in the physical space */
59   PetscObjectState matstate;        /* these are needed just in the saddle point case */
60   PetscObjectState matnnzstate;     /* where we are going to use MatZeroRows on pmat */
61   PetscBool        statechanged;
62   PetscBool        check;
63 } KSP_FETIDP;
64 
65 static PetscErrorCode KSPFETIDPSetPressureOperator_FETIDP(KSP ksp, Mat P) {
66   KSP_FETIDP *fetidp = (KSP_FETIDP *)ksp->data;
67 
68   PetscFunctionBegin;
69   if (P) fetidp->saddlepoint = PETSC_TRUE;
70   PetscCall(PetscObjectCompose((PetscObject)fetidp->innerbddc, "__KSPFETIDP_PPmat", (PetscObject)P));
71   PetscFunctionReturn(0);
72 }
73 
74 /*@
75  KSPFETIDPSetPressureOperator - Sets the operator used to setup the pressure preconditioner for saddle point FETI-DP.
76 
77    Collective on ksp
78 
79    Input Parameters:
80 +  ksp - the FETI-DP Krylov solver
81 -  P - the linear operator to be preconditioned, usually the mass matrix.
82 
83    Level: advanced
84 
85    Notes:
86     The operator can be either passed in a) monolithic global ordering, b) pressure-only global ordering
87           or c) interface pressure ordering (if -ksp_fetidp_pressure_all false).
88           In cases b) and c), the pressure ordering of dofs needs to satisfy
89              pid_1 < pid_2  iff  gid_1 < gid_2
90           where pid_1 and pid_2 are two different pressure dof numbers and gid_1 and gid_2 the corresponding
91           id in the monolithic global ordering.
92 
93 .seealso: `MATIS`, `PCBDDC`, `KSPFETIDPGetInnerBDDC`, `KSPFETIDPGetInnerKSP`, `KSPSetOperators`
94 @*/
95 PetscErrorCode KSPFETIDPSetPressureOperator(KSP ksp, Mat P) {
96   PetscFunctionBegin;
97   PetscValidHeaderSpecific(ksp, KSP_CLASSID, 1);
98   if (P) PetscValidHeaderSpecific(P, MAT_CLASSID, 2);
99   PetscTryMethod(ksp, "KSPFETIDPSetPressureOperator_C", (KSP, Mat), (ksp, P));
100   PetscFunctionReturn(0);
101 }
102 
103 static PetscErrorCode KSPFETIDPGetInnerKSP_FETIDP(KSP ksp, KSP *innerksp) {
104   KSP_FETIDP *fetidp = (KSP_FETIDP *)ksp->data;
105 
106   PetscFunctionBegin;
107   *innerksp = fetidp->innerksp;
108   PetscFunctionReturn(0);
109 }
110 
111 /*@
112  KSPFETIDPGetInnerKSP - Gets the KSP object for the Lagrange multipliers
113 
114    Input Parameters:
115 +  ksp - the FETI-DP KSP
116 -  innerksp - the KSP for the multipliers
117 
118    Level: advanced
119 
120    Notes:
121 
122 .seealso: `MATIS`, `PCBDDC`, `KSPFETIDPSetInnerBDDC`, `KSPFETIDPGetInnerBDDC`
123 @*/
124 PetscErrorCode KSPFETIDPGetInnerKSP(KSP ksp, KSP *innerksp) {
125   PetscFunctionBegin;
126   PetscValidHeaderSpecific(ksp, KSP_CLASSID, 1);
127   PetscValidPointer(innerksp, 2);
128   PetscUseMethod(ksp, "KSPFETIDPGetInnerKSP_C", (KSP, KSP *), (ksp, innerksp));
129   PetscFunctionReturn(0);
130 }
131 
132 static PetscErrorCode KSPFETIDPGetInnerBDDC_FETIDP(KSP ksp, PC *pc) {
133   KSP_FETIDP *fetidp = (KSP_FETIDP *)ksp->data;
134 
135   PetscFunctionBegin;
136   *pc = fetidp->innerbddc;
137   PetscFunctionReturn(0);
138 }
139 
140 /*@
141  KSPFETIDPGetInnerBDDC - Gets the BDDC preconditioner used to setup the FETI-DP matrix for the Lagrange multipliers
142 
143    Input Parameters:
144 +  ksp - the FETI-DP Krylov solver
145 -  pc - the BDDC preconditioner
146 
147    Level: advanced
148 
149    Notes:
150 
151 .seealso: `MATIS`, `PCBDDC`, `KSPFETIDPSetInnerBDDC`, `KSPFETIDPGetInnerKSP`
152 @*/
153 PetscErrorCode KSPFETIDPGetInnerBDDC(KSP ksp, PC *pc) {
154   PetscFunctionBegin;
155   PetscValidHeaderSpecific(ksp, KSP_CLASSID, 1);
156   PetscValidPointer(pc, 2);
157   PetscUseMethod(ksp, "KSPFETIDPGetInnerBDDC_C", (KSP, PC *), (ksp, pc));
158   PetscFunctionReturn(0);
159 }
160 
161 static PetscErrorCode KSPFETIDPSetInnerBDDC_FETIDP(KSP ksp, PC pc) {
162   KSP_FETIDP *fetidp = (KSP_FETIDP *)ksp->data;
163 
164   PetscFunctionBegin;
165   PetscCall(PetscObjectReference((PetscObject)pc));
166   PetscCall(PCDestroy(&fetidp->innerbddc));
167   fetidp->innerbddc = pc;
168   fetidp->userbddc  = PETSC_TRUE;
169   PetscFunctionReturn(0);
170 }
171 
172 /*@
173  KSPFETIDPSetInnerBDDC - Sets the BDDC preconditioner used to setup the FETI-DP matrix for the Lagrange multipliers
174 
175    Collective on ksp
176 
177    Input Parameters:
178 +  ksp - the FETI-DP Krylov solver
179 -  pc - the BDDC preconditioner
180 
181    Level: advanced
182 
183    Notes:
184 
185 .seealso: `MATIS`, `PCBDDC`, `KSPFETIDPGetInnerBDDC`, `KSPFETIDPGetInnerKSP`
186 @*/
187 PetscErrorCode KSPFETIDPSetInnerBDDC(KSP ksp, PC pc) {
188   PetscBool isbddc;
189 
190   PetscFunctionBegin;
191   PetscValidHeaderSpecific(ksp, KSP_CLASSID, 1);
192   PetscValidHeaderSpecific(pc, PC_CLASSID, 2);
193   PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCBDDC, &isbddc));
194   PetscCheck(isbddc, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_WRONG, "KSPFETIDPSetInnerBDDC need a PCBDDC preconditioner");
195   PetscTryMethod(ksp, "KSPFETIDPSetInnerBDDC_C", (KSP, PC), (ksp, pc));
196   PetscFunctionReturn(0);
197 }
198 
199 static PetscErrorCode KSPBuildSolution_FETIDP(KSP ksp, Vec v, Vec *V) {
200   KSP_FETIDP *fetidp = (KSP_FETIDP *)ksp->data;
201   Mat         F;
202   Vec         Xl;
203 
204   PetscFunctionBegin;
205   PetscCall(KSPGetOperators(fetidp->innerksp, &F, NULL));
206   PetscCall(KSPBuildSolution(fetidp->innerksp, NULL, &Xl));
207   if (v) {
208     PetscCall(PCBDDCMatFETIDPGetSolution(F, Xl, v));
209     *V = v;
210   } else {
211     PetscCall(PCBDDCMatFETIDPGetSolution(F, Xl, *V));
212   }
213   PetscFunctionReturn(0);
214 }
215 
216 static PetscErrorCode KSPMonitor_FETIDP(KSP ksp, PetscInt it, PetscReal rnorm, void *ctx) {
217   KSP_FETIDPMon *monctx = (KSP_FETIDPMon *)ctx;
218 
219   PetscFunctionBegin;
220   PetscCall(KSPMonitor(monctx->parentksp, it, rnorm));
221   PetscFunctionReturn(0);
222 }
223 
224 static PetscErrorCode KSPComputeEigenvalues_FETIDP(KSP ksp, PetscInt nmax, PetscReal *r, PetscReal *c, PetscInt *neig) {
225   KSP_FETIDP *fetidp = (KSP_FETIDP *)ksp->data;
226 
227   PetscFunctionBegin;
228   PetscCall(KSPComputeEigenvalues(fetidp->innerksp, nmax, r, c, neig));
229   PetscFunctionReturn(0);
230 }
231 
232 static PetscErrorCode KSPComputeExtremeSingularValues_FETIDP(KSP ksp, PetscReal *emax, PetscReal *emin) {
233   KSP_FETIDP *fetidp = (KSP_FETIDP *)ksp->data;
234 
235   PetscFunctionBegin;
236   PetscCall(KSPComputeExtremeSingularValues(fetidp->innerksp, emax, emin));
237   PetscFunctionReturn(0);
238 }
239 
240 static PetscErrorCode KSPFETIDPCheckOperators(KSP ksp, PetscViewer viewer) {
241   KSP_FETIDP     *fetidp = (KSP_FETIDP *)ksp->data;
242   PC_BDDC        *pcbddc = (PC_BDDC *)fetidp->innerbddc->data;
243   PC_IS          *pcis   = (PC_IS *)fetidp->innerbddc->data;
244   Mat_IS         *matis  = (Mat_IS *)fetidp->innerbddc->pmat->data;
245   Mat             F;
246   FETIDPMat_ctx   fetidpmat_ctx;
247   Vec             test_vec, test_vec_p = NULL, fetidp_global;
248   IS              dirdofs, isvert;
249   MPI_Comm        comm = PetscObjectComm((PetscObject)ksp);
250   PetscScalar     sval, *array;
251   PetscReal       val, rval;
252   const PetscInt *vertex_indices;
253   PetscInt        i, n_vertices;
254   PetscBool       isascii;
255 
256   PetscFunctionBegin;
257   PetscCheckSameComm(ksp, 1, viewer, 2);
258   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
259   PetscCheck(isascii, comm, PETSC_ERR_SUP, "Unsupported viewer");
260   PetscCall(PetscViewerASCIIPrintf(viewer, "----------FETI-DP MAT  --------------\n"));
261   PetscCall(PetscViewerASCIIAddTab(viewer, 2));
262   PetscCall(KSPGetOperators(fetidp->innerksp, &F, NULL));
263   PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_INFO));
264   PetscCall(MatView(F, viewer));
265   PetscCall(PetscViewerPopFormat(viewer));
266   PetscCall(PetscViewerASCIISubtractTab(viewer, 2));
267   PetscCall(MatShellGetContext(F, &fetidpmat_ctx));
268   PetscCall(PetscViewerASCIIPrintf(viewer, "----------FETI-DP TESTS--------------\n"));
269   PetscCall(PetscViewerASCIIPrintf(viewer, "All tests should return zero!\n"));
270   PetscCall(PetscViewerASCIIPrintf(viewer, "FETIDP MAT context in the "));
271   if (fetidp->fully_redundant) {
272     PetscCall(PetscViewerASCIIPrintf(viewer, "fully redundant case for lagrange multipliers.\n"));
273   } else {
274     PetscCall(PetscViewerASCIIPrintf(viewer, "Non-fully redundant case for lagrange multiplier.\n"));
275   }
276   PetscCall(PetscViewerFlush(viewer));
277 
278   /* Get Vertices used to define the BDDC */
279   PetscCall(PCBDDCGraphGetCandidatesIS(pcbddc->mat_graph, NULL, NULL, NULL, NULL, &isvert));
280   PetscCall(ISGetLocalSize(isvert, &n_vertices));
281   PetscCall(ISGetIndices(isvert, &vertex_indices));
282 
283   /******************************************************************/
284   /* TEST A/B: Test numbering of global fetidp dofs                 */
285   /******************************************************************/
286   PetscCall(MatCreateVecs(F, &fetidp_global, NULL));
287   PetscCall(VecDuplicate(fetidpmat_ctx->lambda_local, &test_vec));
288   PetscCall(VecSet(fetidp_global, 1.0));
289   PetscCall(VecSet(test_vec, 1.));
290   PetscCall(VecScatterBegin(fetidpmat_ctx->l2g_lambda, fetidp_global, fetidpmat_ctx->lambda_local, INSERT_VALUES, SCATTER_REVERSE));
291   PetscCall(VecScatterEnd(fetidpmat_ctx->l2g_lambda, fetidp_global, fetidpmat_ctx->lambda_local, INSERT_VALUES, SCATTER_REVERSE));
292   if (fetidpmat_ctx->l2g_p) {
293     PetscCall(VecDuplicate(fetidpmat_ctx->vP, &test_vec_p));
294     PetscCall(VecSet(test_vec_p, 1.));
295     PetscCall(VecScatterBegin(fetidpmat_ctx->l2g_p, fetidp_global, fetidpmat_ctx->vP, INSERT_VALUES, SCATTER_REVERSE));
296     PetscCall(VecScatterEnd(fetidpmat_ctx->l2g_p, fetidp_global, fetidpmat_ctx->vP, INSERT_VALUES, SCATTER_REVERSE));
297   }
298   PetscCall(VecAXPY(test_vec, -1.0, fetidpmat_ctx->lambda_local));
299   PetscCall(VecNorm(test_vec, NORM_INFINITY, &val));
300   PetscCall(VecDestroy(&test_vec));
301   PetscCallMPI(MPI_Reduce(&val, &rval, 1, MPIU_REAL, MPIU_MAX, 0, comm));
302   PetscCall(PetscViewerASCIIPrintf(viewer, "A: CHECK glob to loc: % 1.14e\n", (double)rval));
303 
304   if (fetidpmat_ctx->l2g_p) {
305     PetscCall(VecAXPY(test_vec_p, -1.0, fetidpmat_ctx->vP));
306     PetscCall(VecNorm(test_vec_p, NORM_INFINITY, &val));
307     PetscCallMPI(MPI_Reduce(&val, &rval, 1, MPIU_REAL, MPIU_MAX, 0, comm));
308     PetscCall(PetscViewerASCIIPrintf(viewer, "A: CHECK glob to loc (p): % 1.14e\n", (double)rval));
309   }
310 
311   if (fetidp->fully_redundant) {
312     PetscCall(VecSet(fetidp_global, 0.0));
313     PetscCall(VecSet(fetidpmat_ctx->lambda_local, 0.5));
314     PetscCall(VecScatterBegin(fetidpmat_ctx->l2g_lambda, fetidpmat_ctx->lambda_local, fetidp_global, ADD_VALUES, SCATTER_FORWARD));
315     PetscCall(VecScatterEnd(fetidpmat_ctx->l2g_lambda, fetidpmat_ctx->lambda_local, fetidp_global, ADD_VALUES, SCATTER_FORWARD));
316     PetscCall(VecSum(fetidp_global, &sval));
317     val = PetscRealPart(sval) - fetidpmat_ctx->n_lambda;
318     PetscCallMPI(MPI_Reduce(&val, &rval, 1, MPIU_REAL, MPIU_MAX, 0, comm));
319     PetscCall(PetscViewerASCIIPrintf(viewer, "B: CHECK loc to glob: % 1.14e\n", (double)rval));
320   }
321 
322   if (fetidpmat_ctx->l2g_p) {
323     PetscCall(VecSet(pcis->vec1_N, 1.0));
324     PetscCall(VecSet(pcis->vec1_global, 0.0));
325     PetscCall(VecScatterBegin(matis->rctx, pcis->vec1_N, pcis->vec1_global, ADD_VALUES, SCATTER_REVERSE));
326     PetscCall(VecScatterEnd(matis->rctx, pcis->vec1_N, pcis->vec1_global, ADD_VALUES, SCATTER_REVERSE));
327 
328     PetscCall(VecSet(fetidp_global, 0.0));
329     PetscCall(VecSet(fetidpmat_ctx->vP, -1.0));
330     PetscCall(VecScatterBegin(fetidpmat_ctx->l2g_p, fetidpmat_ctx->vP, fetidp_global, ADD_VALUES, SCATTER_FORWARD));
331     PetscCall(VecScatterEnd(fetidpmat_ctx->l2g_p, fetidpmat_ctx->vP, fetidp_global, ADD_VALUES, SCATTER_FORWARD));
332     PetscCall(VecScatterBegin(fetidpmat_ctx->g2g_p, fetidp_global, pcis->vec1_global, ADD_VALUES, SCATTER_REVERSE));
333     PetscCall(VecScatterEnd(fetidpmat_ctx->g2g_p, fetidp_global, pcis->vec1_global, ADD_VALUES, SCATTER_REVERSE));
334     PetscCall(VecScatterBegin(fetidpmat_ctx->g2g_p, pcis->vec1_global, fetidp_global, INSERT_VALUES, SCATTER_FORWARD));
335     PetscCall(VecScatterEnd(fetidpmat_ctx->g2g_p, pcis->vec1_global, fetidp_global, INSERT_VALUES, SCATTER_FORWARD));
336     PetscCall(VecSum(fetidp_global, &sval));
337     val = PetscRealPart(sval);
338     PetscCallMPI(MPI_Reduce(&val, &rval, 1, MPIU_REAL, MPIU_MAX, 0, comm));
339     PetscCall(PetscViewerASCIIPrintf(viewer, "B: CHECK loc to glob (p): % 1.14e\n", (double)rval));
340   }
341 
342   /******************************************************************/
343   /* TEST C: It should hold B_delta*w=0, w\in\widehat{W}            */
344   /* This is the meaning of the B matrix                            */
345   /******************************************************************/
346 
347   PetscCall(VecSetRandom(pcis->vec1_N, NULL));
348   PetscCall(VecSet(pcis->vec1_global, 0.0));
349   PetscCall(VecScatterBegin(matis->rctx, pcis->vec1_N, pcis->vec1_global, ADD_VALUES, SCATTER_REVERSE));
350   PetscCall(VecScatterEnd(matis->rctx, pcis->vec1_N, pcis->vec1_global, ADD_VALUES, SCATTER_REVERSE));
351   PetscCall(VecScatterBegin(matis->rctx, pcis->vec1_global, pcis->vec1_N, INSERT_VALUES, SCATTER_FORWARD));
352   PetscCall(VecScatterEnd(matis->rctx, pcis->vec1_global, pcis->vec1_N, INSERT_VALUES, SCATTER_FORWARD));
353   PetscCall(VecScatterBegin(pcis->N_to_B, pcis->vec1_N, pcis->vec1_B, INSERT_VALUES, SCATTER_FORWARD));
354   PetscCall(VecScatterEnd(pcis->N_to_B, pcis->vec1_N, pcis->vec1_B, INSERT_VALUES, SCATTER_FORWARD));
355   /* Action of B_delta */
356   PetscCall(MatMult(fetidpmat_ctx->B_delta, pcis->vec1_B, fetidpmat_ctx->lambda_local));
357   PetscCall(VecSet(fetidp_global, 0.0));
358   PetscCall(VecScatterBegin(fetidpmat_ctx->l2g_lambda, fetidpmat_ctx->lambda_local, fetidp_global, ADD_VALUES, SCATTER_FORWARD));
359   PetscCall(VecScatterEnd(fetidpmat_ctx->l2g_lambda, fetidpmat_ctx->lambda_local, fetidp_global, ADD_VALUES, SCATTER_FORWARD));
360   PetscCall(VecNorm(fetidp_global, NORM_INFINITY, &val));
361   PetscCall(PetscViewerASCIIPrintf(viewer, "C: CHECK infty norm of B_delta*w (w continuous): % 1.14e\n", (double)val));
362 
363   /******************************************************************/
364   /* TEST D: It should hold E_Dw = w - P_Dw w\in\widetilde{W}       */
365   /* E_D = R_D^TR                                                   */
366   /* P_D = B_{D,delta}^T B_{delta}                                  */
367   /* eq.44 Mandel Tezaur and Dohrmann 2005                          */
368   /******************************************************************/
369 
370   /* compute a random vector in \widetilde{W} */
371   PetscCall(VecSetRandom(pcis->vec1_N, NULL));
372   /* set zero at vertices and essential dofs */
373   PetscCall(VecGetArray(pcis->vec1_N, &array));
374   for (i = 0; i < n_vertices; i++) array[vertex_indices[i]] = 0.0;
375   PetscCall(PCBDDCGraphGetDirichletDofs(pcbddc->mat_graph, &dirdofs));
376   if (dirdofs) {
377     const PetscInt *idxs;
378     PetscInt        ndir;
379 
380     PetscCall(ISGetLocalSize(dirdofs, &ndir));
381     PetscCall(ISGetIndices(dirdofs, &idxs));
382     for (i = 0; i < ndir; i++) array[idxs[i]] = 0.0;
383     PetscCall(ISRestoreIndices(dirdofs, &idxs));
384   }
385   PetscCall(VecRestoreArray(pcis->vec1_N, &array));
386   /* store w for final comparison */
387   PetscCall(VecDuplicate(pcis->vec1_B, &test_vec));
388   PetscCall(VecScatterBegin(pcis->N_to_B, pcis->vec1_N, test_vec, INSERT_VALUES, SCATTER_FORWARD));
389   PetscCall(VecScatterEnd(pcis->N_to_B, pcis->vec1_N, test_vec, INSERT_VALUES, SCATTER_FORWARD));
390 
391   /* Jump operator P_D : results stored in pcis->vec1_B */
392   /* Action of B_delta */
393   PetscCall(MatMult(fetidpmat_ctx->B_delta, test_vec, fetidpmat_ctx->lambda_local));
394   PetscCall(VecSet(fetidp_global, 0.0));
395   PetscCall(VecScatterBegin(fetidpmat_ctx->l2g_lambda, fetidpmat_ctx->lambda_local, fetidp_global, ADD_VALUES, SCATTER_FORWARD));
396   PetscCall(VecScatterEnd(fetidpmat_ctx->l2g_lambda, fetidpmat_ctx->lambda_local, fetidp_global, ADD_VALUES, SCATTER_FORWARD));
397   /* Action of B_Ddelta^T */
398   PetscCall(VecScatterBegin(fetidpmat_ctx->l2g_lambda, fetidp_global, fetidpmat_ctx->lambda_local, INSERT_VALUES, SCATTER_REVERSE));
399   PetscCall(VecScatterEnd(fetidpmat_ctx->l2g_lambda, fetidp_global, fetidpmat_ctx->lambda_local, INSERT_VALUES, SCATTER_REVERSE));
400   PetscCall(MatMultTranspose(fetidpmat_ctx->B_Ddelta, fetidpmat_ctx->lambda_local, pcis->vec1_B));
401 
402   /* Average operator E_D : results stored in pcis->vec2_B */
403   PetscCall(PCBDDCScalingExtension(fetidpmat_ctx->pc, test_vec, pcis->vec1_global));
404   PetscCall(VecScatterBegin(pcis->global_to_B, pcis->vec1_global, pcis->vec2_B, INSERT_VALUES, SCATTER_FORWARD));
405   PetscCall(VecScatterEnd(pcis->global_to_B, pcis->vec1_global, pcis->vec2_B, INSERT_VALUES, SCATTER_FORWARD));
406 
407   /* test E_D=I-P_D */
408   PetscCall(VecAXPY(pcis->vec1_B, 1.0, pcis->vec2_B));
409   PetscCall(VecAXPY(pcis->vec1_B, -1.0, test_vec));
410   PetscCall(VecNorm(pcis->vec1_B, NORM_INFINITY, &val));
411   PetscCall(VecDestroy(&test_vec));
412   PetscCallMPI(MPI_Reduce(&val, &rval, 1, MPIU_REAL, MPIU_MAX, 0, comm));
413   PetscCall(PetscViewerASCIIPrintf(viewer, "%d: CHECK infty norm of E_D + P_D - I: %1.14e\n", PetscGlobalRank, (double)val));
414 
415   /******************************************************************/
416   /* TEST E: It should hold R_D^TP_Dw=0 w\in\widetilde{W}           */
417   /* eq.48 Mandel Tezaur and Dohrmann 2005                          */
418   /******************************************************************/
419 
420   PetscCall(VecSetRandom(pcis->vec1_N, NULL));
421   /* set zero at vertices and essential dofs */
422   PetscCall(VecGetArray(pcis->vec1_N, &array));
423   for (i = 0; i < n_vertices; i++) array[vertex_indices[i]] = 0.0;
424   if (dirdofs) {
425     const PetscInt *idxs;
426     PetscInt        ndir;
427 
428     PetscCall(ISGetLocalSize(dirdofs, &ndir));
429     PetscCall(ISGetIndices(dirdofs, &idxs));
430     for (i = 0; i < ndir; i++) array[idxs[i]] = 0.0;
431     PetscCall(ISRestoreIndices(dirdofs, &idxs));
432   }
433   PetscCall(VecRestoreArray(pcis->vec1_N, &array));
434 
435   /* Jump operator P_D : results stored in pcis->vec1_B */
436 
437   PetscCall(VecScatterBegin(pcis->N_to_B, pcis->vec1_N, pcis->vec1_B, INSERT_VALUES, SCATTER_FORWARD));
438   PetscCall(VecScatterEnd(pcis->N_to_B, pcis->vec1_N, pcis->vec1_B, INSERT_VALUES, SCATTER_FORWARD));
439   /* Action of B_delta */
440   PetscCall(MatMult(fetidpmat_ctx->B_delta, pcis->vec1_B, fetidpmat_ctx->lambda_local));
441   PetscCall(VecSet(fetidp_global, 0.0));
442   PetscCall(VecScatterBegin(fetidpmat_ctx->l2g_lambda, fetidpmat_ctx->lambda_local, fetidp_global, ADD_VALUES, SCATTER_FORWARD));
443   PetscCall(VecScatterEnd(fetidpmat_ctx->l2g_lambda, fetidpmat_ctx->lambda_local, fetidp_global, ADD_VALUES, SCATTER_FORWARD));
444   /* Action of B_Ddelta^T */
445   PetscCall(VecScatterBegin(fetidpmat_ctx->l2g_lambda, fetidp_global, fetidpmat_ctx->lambda_local, INSERT_VALUES, SCATTER_REVERSE));
446   PetscCall(VecScatterEnd(fetidpmat_ctx->l2g_lambda, fetidp_global, fetidpmat_ctx->lambda_local, INSERT_VALUES, SCATTER_REVERSE));
447   PetscCall(MatMultTranspose(fetidpmat_ctx->B_Ddelta, fetidpmat_ctx->lambda_local, pcis->vec1_B));
448   /* scaling */
449   PetscCall(PCBDDCScalingExtension(fetidpmat_ctx->pc, pcis->vec1_B, pcis->vec1_global));
450   PetscCall(VecNorm(pcis->vec1_global, NORM_INFINITY, &val));
451   PetscCall(PetscViewerASCIIPrintf(viewer, "E: CHECK infty norm of R^T_D P_D: % 1.14e\n", (double)val));
452 
453   if (!fetidp->fully_redundant) {
454     /******************************************************************/
455     /* TEST F: It should holds B_{delta}B^T_{D,delta}=I               */
456     /* Corollary thm 14 Mandel Tezaur and Dohrmann 2005               */
457     /******************************************************************/
458     PetscCall(VecDuplicate(fetidp_global, &test_vec));
459     PetscCall(VecSetRandom(fetidp_global, NULL));
460     if (fetidpmat_ctx->l2g_p) {
461       PetscCall(VecSet(fetidpmat_ctx->vP, 0.));
462       PetscCall(VecScatterBegin(fetidpmat_ctx->l2g_p, fetidpmat_ctx->vP, fetidp_global, INSERT_VALUES, SCATTER_FORWARD));
463       PetscCall(VecScatterEnd(fetidpmat_ctx->l2g_p, fetidpmat_ctx->vP, fetidp_global, INSERT_VALUES, SCATTER_FORWARD));
464     }
465     /* Action of B_Ddelta^T */
466     PetscCall(VecScatterBegin(fetidpmat_ctx->l2g_lambda, fetidp_global, fetidpmat_ctx->lambda_local, INSERT_VALUES, SCATTER_REVERSE));
467     PetscCall(VecScatterEnd(fetidpmat_ctx->l2g_lambda, fetidp_global, fetidpmat_ctx->lambda_local, INSERT_VALUES, SCATTER_REVERSE));
468     PetscCall(MatMultTranspose(fetidpmat_ctx->B_Ddelta, fetidpmat_ctx->lambda_local, pcis->vec1_B));
469     /* Action of B_delta */
470     PetscCall(MatMult(fetidpmat_ctx->B_delta, pcis->vec1_B, fetidpmat_ctx->lambda_local));
471     PetscCall(VecSet(test_vec, 0.0));
472     PetscCall(VecScatterBegin(fetidpmat_ctx->l2g_lambda, fetidpmat_ctx->lambda_local, test_vec, ADD_VALUES, SCATTER_FORWARD));
473     PetscCall(VecScatterEnd(fetidpmat_ctx->l2g_lambda, fetidpmat_ctx->lambda_local, test_vec, ADD_VALUES, SCATTER_FORWARD));
474     PetscCall(VecAXPY(fetidp_global, -1., test_vec));
475     PetscCall(VecNorm(fetidp_global, NORM_INFINITY, &val));
476     PetscCall(PetscViewerASCIIPrintf(viewer, "E: CHECK infty norm of P^T_D - I: % 1.14e\n", (double)val));
477     PetscCall(VecDestroy(&test_vec));
478   }
479   PetscCall(PetscViewerASCIIPrintf(viewer, "-------------------------------------\n"));
480   PetscCall(PetscViewerFlush(viewer));
481   PetscCall(VecDestroy(&test_vec_p));
482   PetscCall(ISDestroy(&dirdofs));
483   PetscCall(VecDestroy(&fetidp_global));
484   PetscCall(ISRestoreIndices(isvert, &vertex_indices));
485   PetscCall(PCBDDCGraphRestoreCandidatesIS(pcbddc->mat_graph, NULL, NULL, NULL, NULL, &isvert));
486   PetscFunctionReturn(0);
487 }
488 
489 static PetscErrorCode KSPFETIDPSetUpOperators(KSP ksp) {
490   KSP_FETIDP      *fetidp = (KSP_FETIDP *)ksp->data;
491   PC_BDDC         *pcbddc = (PC_BDDC *)fetidp->innerbddc->data;
492   Mat              A, Ap;
493   PetscInt         fid = -1;
494   PetscMPIInt      size;
495   PetscBool        ismatis, pisz = PETSC_FALSE, allp = PETSC_FALSE, schp = PETSC_FALSE;
496   PetscBool        flip = PETSC_FALSE; /* Usually, Stokes is written (B = -\int_\Omega \nabla \cdot u q)
497                            | A B'| | v | = | f |
498                            | B 0 | | p | = | g |
499                             If -ksp_fetidp_saddlepoint_flip is true, the code assumes it is written as
500                            | A B'| | v | = | f |
501                            |-B 0 | | p | = |-g |
502                          */
503   PetscObjectState matstate, matnnzstate;
504 
505   PetscFunctionBegin;
506   PetscOptionsBegin(PetscObjectComm((PetscObject)ksp), ((PetscObject)ksp)->prefix, "FETI-DP options", "PC");
507   PetscCall(PetscOptionsInt("-ksp_fetidp_pressure_field", "Field id for pressures for saddle-point problems", NULL, fid, &fid, NULL));
508   PetscCall(PetscOptionsBool("-ksp_fetidp_pressure_all", "Use the whole pressure set instead of just that at the interface", NULL, allp, &allp, NULL));
509   PetscCall(PetscOptionsBool("-ksp_fetidp_saddlepoint_flip", "Flip the sign of the pressure-velocity (lower-left) block", NULL, flip, &flip, NULL));
510   PetscCall(PetscOptionsBool("-ksp_fetidp_pressure_schur", "Use a BDDC solver for pressure", NULL, schp, &schp, NULL));
511   PetscOptionsEnd();
512 
513   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)ksp), &size));
514   fetidp->saddlepoint = (fid >= 0 ? PETSC_TRUE : fetidp->saddlepoint);
515   if (size == 1) fetidp->saddlepoint = PETSC_FALSE;
516 
517   PetscCall(KSPGetOperators(ksp, &A, &Ap));
518   PetscCall(PetscObjectTypeCompare((PetscObject)A, MATIS, &ismatis));
519   PetscCheck(ismatis, PetscObjectComm((PetscObject)ksp), PETSC_ERR_USER, "Amat should be of type MATIS");
520 
521   /* Quiet return if the matrix states are unchanged.
522      Needed only for the saddle point case since it uses MatZeroRows
523      on a matrix that may not have changed */
524   PetscCall(PetscObjectStateGet((PetscObject)A, &matstate));
525   PetscCall(MatGetNonzeroState(A, &matnnzstate));
526   if (matstate == fetidp->matstate && matnnzstate == fetidp->matnnzstate) PetscFunctionReturn(0);
527   fetidp->matstate     = matstate;
528   fetidp->matnnzstate  = matnnzstate;
529   fetidp->statechanged = fetidp->saddlepoint;
530 
531   /* see if we have some fields attached */
532   if (!pcbddc->n_ISForDofsLocal && !pcbddc->n_ISForDofs) {
533     DM             dm;
534     PetscContainer c;
535 
536     PetscCall(KSPGetDM(ksp, &dm));
537     PetscCall(PetscObjectQuery((PetscObject)A, "_convert_nest_lfields", (PetscObject *)&c));
538     if (dm) {
539       IS      *fields;
540       PetscInt nf, i;
541 
542       PetscCall(DMCreateFieldDecomposition(dm, &nf, NULL, &fields, NULL));
543       PetscCall(PCBDDCSetDofsSplitting(fetidp->innerbddc, nf, fields));
544       for (i = 0; i < nf; i++) { PetscCall(ISDestroy(&fields[i])); }
545       PetscCall(PetscFree(fields));
546     } else if (c) {
547       MatISLocalFields lf;
548 
549       PetscCall(PetscContainerGetPointer(c, (void **)&lf));
550       PetscCall(PCBDDCSetDofsSplittingLocal(fetidp->innerbddc, lf->nr, lf->rf));
551     }
552   }
553 
554   if (!fetidp->saddlepoint) {
555     PetscCall(PCSetOperators(fetidp->innerbddc, A, A));
556   } else {
557     Mat          nA, lA, PPmat;
558     MatNullSpace nnsp;
559     IS           pP;
560     PetscInt     totP;
561 
562     PetscCall(MatISGetLocalMat(A, &lA));
563     PetscCall(PetscObjectCompose((PetscObject)fetidp->innerbddc, "__KSPFETIDP_lA", (PetscObject)lA));
564 
565     pP = fetidp->pP;
566     if (!pP) { /* first time, need to compute pressure dofs */
567       PC_IS                 *pcis  = (PC_IS *)fetidp->innerbddc->data;
568       Mat_IS                *matis = (Mat_IS *)(A->data);
569       ISLocalToGlobalMapping l2g;
570       IS                     lP = NULL, II, pII, lPall, Pall, is1, is2;
571       const PetscInt        *idxs;
572       PetscInt               nl, ni, *widxs;
573       PetscInt               i, j, n_neigh, *neigh, *n_shared, **shared, *count;
574       PetscInt               rst, ren, n;
575       PetscBool              ploc;
576 
577       PetscCall(MatGetLocalSize(A, &nl, NULL));
578       PetscCall(MatGetOwnershipRange(A, &rst, &ren));
579       PetscCall(MatGetLocalSize(lA, &n, NULL));
580       PetscCall(MatISGetLocalToGlobalMapping(A, &l2g, NULL));
581 
582       if (!pcis->is_I_local) { /* need to compute interior dofs */
583         PetscCall(PetscCalloc1(n, &count));
584         PetscCall(ISLocalToGlobalMappingGetInfo(l2g, &n_neigh, &neigh, &n_shared, &shared));
585         for (i = 1; i < n_neigh; i++)
586           for (j = 0; j < n_shared[i]; j++) count[shared[i][j]] += 1;
587         for (i = 0, j = 0; i < n; i++)
588           if (!count[i]) count[j++] = i;
589         PetscCall(ISLocalToGlobalMappingRestoreInfo(l2g, &n_neigh, &neigh, &n_shared, &shared));
590         PetscCall(ISCreateGeneral(PETSC_COMM_SELF, j, count, PETSC_OWN_POINTER, &II));
591       } else {
592         PetscCall(PetscObjectReference((PetscObject)pcis->is_I_local));
593         II = pcis->is_I_local;
594       }
595 
596       /* interior dofs in layout */
597       PetscCall(PetscArrayzero(matis->sf_leafdata, n));
598       PetscCall(PetscArrayzero(matis->sf_rootdata, nl));
599       PetscCall(ISGetLocalSize(II, &ni));
600       PetscCall(ISGetIndices(II, &idxs));
601       for (i = 0; i < ni; i++) matis->sf_leafdata[idxs[i]] = 1;
602       PetscCall(ISRestoreIndices(II, &idxs));
603       PetscCall(PetscSFReduceBegin(matis->sf, MPIU_INT, matis->sf_leafdata, matis->sf_rootdata, MPI_REPLACE));
604       PetscCall(PetscSFReduceEnd(matis->sf, MPIU_INT, matis->sf_leafdata, matis->sf_rootdata, MPI_REPLACE));
605       PetscCall(PetscMalloc1(PetscMax(nl, n), &widxs));
606       for (i = 0, ni = 0; i < nl; i++)
607         if (matis->sf_rootdata[i]) widxs[ni++] = i + rst;
608       PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)ksp), ni, widxs, PETSC_COPY_VALUES, &pII));
609 
610       /* pressure dofs */
611       Pall  = NULL;
612       lPall = NULL;
613       ploc  = PETSC_FALSE;
614       if (fid < 0) { /* zero pressure block */
615         PetscInt np;
616 
617         PetscCall(MatFindZeroDiagonals(A, &Pall));
618         PetscCall(ISGetSize(Pall, &np));
619         if (!np) { /* zero-block not found, defaults to last field (if set) */
620           fid = pcbddc->n_ISForDofsLocal ? pcbddc->n_ISForDofsLocal - 1 : pcbddc->n_ISForDofs - 1;
621           PetscCall(ISDestroy(&Pall));
622         } else if (!pcbddc->n_ISForDofsLocal && !pcbddc->n_ISForDofs) {
623           PetscCall(PCBDDCSetDofsSplitting(fetidp->innerbddc, 1, &Pall));
624         }
625       }
626       if (!Pall) { /* look for registered fields */
627         if (pcbddc->n_ISForDofsLocal) {
628           PetscInt np;
629 
630           PetscCheck(fid >= 0 && fid < pcbddc->n_ISForDofsLocal, PetscObjectComm((PetscObject)ksp), PETSC_ERR_USER, "Invalid field id for pressure %" PetscInt_FMT ", max %" PetscInt_FMT, fid, pcbddc->n_ISForDofsLocal);
631           /* need a sequential IS */
632           PetscCall(ISGetLocalSize(pcbddc->ISForDofsLocal[fid], &np));
633           PetscCall(ISGetIndices(pcbddc->ISForDofsLocal[fid], &idxs));
634           PetscCall(ISCreateGeneral(PETSC_COMM_SELF, np, idxs, PETSC_COPY_VALUES, &lPall));
635           PetscCall(ISRestoreIndices(pcbddc->ISForDofsLocal[fid], &idxs));
636           ploc = PETSC_TRUE;
637         } else if (pcbddc->n_ISForDofs) {
638           PetscCheck(fid >= 0 && fid < pcbddc->n_ISForDofs, PetscObjectComm((PetscObject)ksp), PETSC_ERR_USER, "Invalid field id for pressure %" PetscInt_FMT ", max %" PetscInt_FMT, fid, pcbddc->n_ISForDofs);
639           PetscCall(PetscObjectReference((PetscObject)pcbddc->ISForDofs[fid]));
640           Pall = pcbddc->ISForDofs[fid];
641         } else SETERRQ(PetscObjectComm((PetscObject)ksp), PETSC_ERR_USER, "Cannot detect pressure field! Use KSPFETIDPGetInnerBDDC() + PCBDDCSetDofsSplitting or PCBDDCSetDofsSplittingLocal");
642       }
643 
644       /* if the user requested the entire pressure,
645          remove the interior pressure dofs from II (or pII) */
646       if (allp) {
647         if (ploc) {
648           IS nII;
649           PetscCall(ISDifference(II, lPall, &nII));
650           PetscCall(ISDestroy(&II));
651           II = nII;
652         } else {
653           IS nII;
654           PetscCall(ISDifference(pII, Pall, &nII));
655           PetscCall(ISDestroy(&pII));
656           pII = nII;
657         }
658       }
659       if (ploc) {
660         PetscCall(ISDifference(lPall, II, &lP));
661         PetscCall(PetscObjectCompose((PetscObject)fetidp->innerbddc, "__KSPFETIDP_lP", (PetscObject)lP));
662       } else {
663         PetscCall(ISDifference(Pall, pII, &pP));
664         PetscCall(PetscObjectCompose((PetscObject)fetidp->innerbddc, "__KSPFETIDP_pP", (PetscObject)pP));
665         /* need all local pressure dofs */
666         PetscCall(PetscArrayzero(matis->sf_leafdata, n));
667         PetscCall(PetscArrayzero(matis->sf_rootdata, nl));
668         PetscCall(ISGetLocalSize(Pall, &ni));
669         PetscCall(ISGetIndices(Pall, &idxs));
670         for (i = 0; i < ni; i++) matis->sf_rootdata[idxs[i] - rst] = 1;
671         PetscCall(ISRestoreIndices(Pall, &idxs));
672         PetscCall(PetscSFBcastBegin(matis->sf, MPIU_INT, matis->sf_rootdata, matis->sf_leafdata, MPI_REPLACE));
673         PetscCall(PetscSFBcastEnd(matis->sf, MPIU_INT, matis->sf_rootdata, matis->sf_leafdata, MPI_REPLACE));
674         for (i = 0, ni = 0; i < n; i++)
675           if (matis->sf_leafdata[i]) widxs[ni++] = i;
676         PetscCall(ISCreateGeneral(PETSC_COMM_SELF, ni, widxs, PETSC_COPY_VALUES, &lPall));
677       }
678 
679       if (!Pall) {
680         PetscCall(PetscArrayzero(matis->sf_leafdata, n));
681         PetscCall(PetscArrayzero(matis->sf_rootdata, nl));
682         PetscCall(ISGetLocalSize(lPall, &ni));
683         PetscCall(ISGetIndices(lPall, &idxs));
684         for (i = 0; i < ni; i++) matis->sf_leafdata[idxs[i]] = 1;
685         PetscCall(ISRestoreIndices(lPall, &idxs));
686         PetscCall(PetscSFReduceBegin(matis->sf, MPIU_INT, matis->sf_leafdata, matis->sf_rootdata, MPI_REPLACE));
687         PetscCall(PetscSFReduceEnd(matis->sf, MPIU_INT, matis->sf_leafdata, matis->sf_rootdata, MPI_REPLACE));
688         for (i = 0, ni = 0; i < nl; i++)
689           if (matis->sf_rootdata[i]) widxs[ni++] = i + rst;
690         PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)ksp), ni, widxs, PETSC_COPY_VALUES, &Pall));
691       }
692       PetscCall(PetscObjectCompose((PetscObject)fetidp->innerbddc, "__KSPFETIDP_aP", (PetscObject)Pall));
693 
694       if (flip) {
695         PetscInt npl;
696         PetscCall(ISGetLocalSize(Pall, &npl));
697         PetscCall(ISGetIndices(Pall, &idxs));
698         PetscCall(MatCreateVecs(A, NULL, &fetidp->rhs_flip));
699         PetscCall(VecSet(fetidp->rhs_flip, 1.));
700         PetscCall(VecSetOption(fetidp->rhs_flip, VEC_IGNORE_OFF_PROC_ENTRIES, PETSC_TRUE));
701         for (i = 0; i < npl; i++) { PetscCall(VecSetValue(fetidp->rhs_flip, idxs[i], -1., INSERT_VALUES)); }
702         PetscCall(VecAssemblyBegin(fetidp->rhs_flip));
703         PetscCall(VecAssemblyEnd(fetidp->rhs_flip));
704         PetscCall(PetscObjectCompose((PetscObject)fetidp->innerbddc, "__KSPFETIDP_flip", (PetscObject)fetidp->rhs_flip));
705         PetscCall(ISRestoreIndices(Pall, &idxs));
706       }
707       PetscCall(ISDestroy(&Pall));
708       PetscCall(ISDestroy(&pII));
709 
710       /* local selected pressures in subdomain-wise and global ordering */
711       PetscCall(PetscArrayzero(matis->sf_leafdata, n));
712       PetscCall(PetscArrayzero(matis->sf_rootdata, nl));
713       if (!ploc) {
714         PetscInt *widxs2;
715 
716         PetscCheck(pP, PetscObjectComm((PetscObject)ksp), PETSC_ERR_PLIB, "Missing parallel pressure IS");
717         PetscCall(ISGetLocalSize(pP, &ni));
718         PetscCall(ISGetIndices(pP, &idxs));
719         for (i = 0; i < ni; i++) matis->sf_rootdata[idxs[i] - rst] = 1;
720         PetscCall(ISRestoreIndices(pP, &idxs));
721         PetscCall(PetscSFBcastBegin(matis->sf, MPIU_INT, matis->sf_rootdata, matis->sf_leafdata, MPI_REPLACE));
722         PetscCall(PetscSFBcastEnd(matis->sf, MPIU_INT, matis->sf_rootdata, matis->sf_leafdata, MPI_REPLACE));
723         for (i = 0, ni = 0; i < n; i++)
724           if (matis->sf_leafdata[i]) widxs[ni++] = i;
725         PetscCall(PetscMalloc1(ni, &widxs2));
726         PetscCall(ISLocalToGlobalMappingApply(l2g, ni, widxs, widxs2));
727         PetscCall(ISCreateGeneral(PETSC_COMM_SELF, ni, widxs, PETSC_COPY_VALUES, &lP));
728         PetscCall(PetscObjectCompose((PetscObject)fetidp->innerbddc, "__KSPFETIDP_lP", (PetscObject)lP));
729         PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)ksp), ni, widxs2, PETSC_OWN_POINTER, &is1));
730         PetscCall(PetscObjectCompose((PetscObject)fetidp->innerbddc, "__KSPFETIDP_gP", (PetscObject)is1));
731         PetscCall(ISDestroy(&is1));
732       } else {
733         PetscCheck(lP, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Missing sequential pressure IS");
734         PetscCall(ISGetLocalSize(lP, &ni));
735         PetscCall(ISGetIndices(lP, &idxs));
736         for (i = 0; i < ni; i++)
737           if (idxs[i] >= 0 && idxs[i] < n) matis->sf_leafdata[idxs[i]] = 1;
738         PetscCall(ISRestoreIndices(lP, &idxs));
739         PetscCall(PetscSFReduceBegin(matis->sf, MPIU_INT, matis->sf_leafdata, matis->sf_rootdata, MPI_REPLACE));
740         PetscCall(ISLocalToGlobalMappingApply(l2g, ni, idxs, widxs));
741         PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)ksp), ni, widxs, PETSC_COPY_VALUES, &is1));
742         PetscCall(PetscObjectCompose((PetscObject)fetidp->innerbddc, "__KSPFETIDP_gP", (PetscObject)is1));
743         PetscCall(ISDestroy(&is1));
744         PetscCall(PetscSFReduceEnd(matis->sf, MPIU_INT, matis->sf_leafdata, matis->sf_rootdata, MPI_REPLACE));
745         for (i = 0, ni = 0; i < nl; i++)
746           if (matis->sf_rootdata[i]) widxs[ni++] = i + rst;
747         PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)ksp), ni, widxs, PETSC_COPY_VALUES, &pP));
748         PetscCall(PetscObjectCompose((PetscObject)fetidp->innerbddc, "__KSPFETIDP_pP", (PetscObject)pP));
749       }
750       PetscCall(PetscFree(widxs));
751 
752       /* If there's any "interior pressure",
753          we may want to use a discrete harmonic solver instead
754          of a Stokes harmonic for the Dirichlet preconditioner
755          Need to extract the interior velocity dofs in interior dofs ordering (iV)
756          and interior pressure dofs in local ordering (iP) */
757       if (!allp) {
758         ISLocalToGlobalMapping l2g_t;
759 
760         PetscCall(ISDifference(lPall, lP, &is1));
761         PetscCall(PetscObjectCompose((PetscObject)fetidp->innerbddc, "__KSPFETIDP_iP", (PetscObject)is1));
762         PetscCall(ISDifference(II, is1, &is2));
763         PetscCall(ISDestroy(&is1));
764         PetscCall(ISLocalToGlobalMappingCreateIS(II, &l2g_t));
765         PetscCall(ISGlobalToLocalMappingApplyIS(l2g_t, IS_GTOLM_DROP, is2, &is1));
766         PetscCall(ISGetLocalSize(is1, &i));
767         PetscCall(ISGetLocalSize(is2, &j));
768         PetscCheck(i == j, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent local sizes %" PetscInt_FMT " and %" PetscInt_FMT " for iV", i, j);
769         PetscCall(PetscObjectCompose((PetscObject)fetidp->innerbddc, "__KSPFETIDP_iV", (PetscObject)is1));
770         PetscCall(ISLocalToGlobalMappingDestroy(&l2g_t));
771         PetscCall(ISDestroy(&is1));
772         PetscCall(ISDestroy(&is2));
773       }
774       PetscCall(ISDestroy(&II));
775 
776       /* exclude selected pressures from the inner BDDC */
777       if (pcbddc->DirichletBoundariesLocal) {
778         IS       list[2], plP, isout;
779         PetscInt np;
780 
781         /* need a parallel IS */
782         PetscCall(ISGetLocalSize(lP, &np));
783         PetscCall(ISGetIndices(lP, &idxs));
784         PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)ksp), np, idxs, PETSC_USE_POINTER, &plP));
785         list[0] = plP;
786         list[1] = pcbddc->DirichletBoundariesLocal;
787         PetscCall(ISConcatenate(PetscObjectComm((PetscObject)ksp), 2, list, &isout));
788         PetscCall(ISSortRemoveDups(isout));
789         PetscCall(ISDestroy(&plP));
790         PetscCall(ISRestoreIndices(lP, &idxs));
791         PetscCall(PCBDDCSetDirichletBoundariesLocal(fetidp->innerbddc, isout));
792         PetscCall(ISDestroy(&isout));
793       } else if (pcbddc->DirichletBoundaries) {
794         IS list[2], isout;
795 
796         list[0] = pP;
797         list[1] = pcbddc->DirichletBoundaries;
798         PetscCall(ISConcatenate(PetscObjectComm((PetscObject)ksp), 2, list, &isout));
799         PetscCall(ISSortRemoveDups(isout));
800         PetscCall(PCBDDCSetDirichletBoundaries(fetidp->innerbddc, isout));
801         PetscCall(ISDestroy(&isout));
802       } else {
803         IS       plP;
804         PetscInt np;
805 
806         /* need a parallel IS */
807         PetscCall(ISGetLocalSize(lP, &np));
808         PetscCall(ISGetIndices(lP, &idxs));
809         PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)ksp), np, idxs, PETSC_COPY_VALUES, &plP));
810         PetscCall(PCBDDCSetDirichletBoundariesLocal(fetidp->innerbddc, plP));
811         PetscCall(ISDestroy(&plP));
812         PetscCall(ISRestoreIndices(lP, &idxs));
813       }
814 
815       /* save CSR information for the pressure BDDC solver (if any) */
816       if (schp) {
817         PetscInt np, nt;
818 
819         PetscCall(MatGetSize(matis->A, &nt, NULL));
820         PetscCall(ISGetLocalSize(lP, &np));
821         if (np) {
822           PetscInt *xadj = pcbddc->mat_graph->xadj;
823           PetscInt *adjn = pcbddc->mat_graph->adjncy;
824           PetscInt  nv   = pcbddc->mat_graph->nvtxs_csr;
825 
826           if (nv && nv == nt) {
827             ISLocalToGlobalMapping pmap;
828             PetscInt              *schp_csr, *schp_xadj, *schp_adjn, p;
829             PetscContainer         c;
830 
831             PetscCall(ISLocalToGlobalMappingCreateIS(lPall, &pmap));
832             PetscCall(ISGetIndices(lPall, &idxs));
833             for (p = 0, nv = 0; p < np; p++) {
834               PetscInt x, n = idxs[p];
835 
836               PetscCall(ISGlobalToLocalMappingApply(pmap, IS_GTOLM_DROP, xadj[n + 1] - xadj[n], adjn + xadj[n], &x, NULL));
837               nv += x;
838             }
839             PetscCall(PetscMalloc1(np + 1 + nv, &schp_csr));
840             schp_xadj = schp_csr;
841             schp_adjn = schp_csr + np + 1;
842             for (p = 0, schp_xadj[0] = 0; p < np; p++) {
843               PetscInt x, n = idxs[p];
844 
845               PetscCall(ISGlobalToLocalMappingApply(pmap, IS_GTOLM_DROP, xadj[n + 1] - xadj[n], adjn + xadj[n], &x, schp_adjn + schp_xadj[p]));
846               schp_xadj[p + 1] = schp_xadj[p] + x;
847             }
848             PetscCall(ISRestoreIndices(lPall, &idxs));
849             PetscCall(ISLocalToGlobalMappingDestroy(&pmap));
850             PetscCall(PetscContainerCreate(PETSC_COMM_SELF, &c));
851             PetscCall(PetscContainerSetPointer(c, schp_csr));
852             PetscCall(PetscContainerSetUserDestroy(c, PetscContainerUserDestroyDefault));
853             PetscCall(PetscObjectCompose((PetscObject)fetidp->innerbddc, "__KSPFETIDP_pCSR", (PetscObject)c));
854             PetscCall(PetscContainerDestroy(&c));
855           }
856         }
857       }
858       PetscCall(ISDestroy(&lPall));
859       PetscCall(ISDestroy(&lP));
860       fetidp->pP = pP;
861     }
862 
863     /* total number of selected pressure dofs */
864     PetscCall(ISGetSize(fetidp->pP, &totP));
865 
866     /* Set operator for inner BDDC */
867     if (totP || fetidp->rhs_flip) {
868       PetscCall(MatDuplicate(A, MAT_COPY_VALUES, &nA));
869     } else {
870       PetscCall(PetscObjectReference((PetscObject)A));
871       nA = A;
872     }
873     if (fetidp->rhs_flip) {
874       PetscCall(MatDiagonalScale(nA, fetidp->rhs_flip, NULL));
875       if (totP) {
876         Mat lA2;
877 
878         PetscCall(MatISGetLocalMat(nA, &lA));
879         PetscCall(MatDuplicate(lA, MAT_COPY_VALUES, &lA2));
880         PetscCall(PetscObjectCompose((PetscObject)fetidp->innerbddc, "__KSPFETIDP_lA", (PetscObject)lA2));
881         PetscCall(MatDestroy(&lA2));
882       }
883     }
884 
885     if (totP) {
886       PetscCall(MatSetOption(nA, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_FALSE));
887       PetscCall(MatZeroRowsColumnsIS(nA, fetidp->pP, 1., NULL, NULL));
888     } else {
889       PetscCall(PetscObjectCompose((PetscObject)fetidp->innerbddc, "__KSPFETIDP_lA", NULL));
890     }
891     PetscCall(MatGetNearNullSpace(Ap, &nnsp));
892     if (!nnsp) { PetscCall(MatGetNullSpace(Ap, &nnsp)); }
893     if (!nnsp) { PetscCall(MatGetNearNullSpace(A, &nnsp)); }
894     if (!nnsp) { PetscCall(MatGetNullSpace(A, &nnsp)); }
895     PetscCall(MatSetNearNullSpace(nA, nnsp));
896     PetscCall(PCSetOperators(fetidp->innerbddc, nA, nA));
897     PetscCall(MatDestroy(&nA));
898 
899     /* non-zero rhs on interior dofs when applying the preconditioner */
900     if (totP) pcbddc->switch_static = PETSC_TRUE;
901 
902     /* if there are no interface pressures, set inner bddc flag for benign saddle point */
903     if (!totP) {
904       pcbddc->benign_saddle_point = PETSC_TRUE;
905       pcbddc->compute_nonetflux   = PETSC_TRUE;
906     }
907 
908     /* Operators for pressure preconditioner */
909     if (totP) {
910       /* Extract pressure block if needed */
911       if (!pisz) {
912         Mat C;
913         IS  nzrows = NULL;
914 
915         PetscCall(MatCreateSubMatrix(A, fetidp->pP, fetidp->pP, MAT_INITIAL_MATRIX, &C));
916         PetscCall(MatFindNonzeroRows(C, &nzrows));
917         if (nzrows) {
918           PetscInt i;
919 
920           PetscCall(ISGetSize(nzrows, &i));
921           PetscCall(ISDestroy(&nzrows));
922           if (!i) pisz = PETSC_TRUE;
923         }
924         if (!pisz) {
925           PetscCall(MatScale(C, -1.)); /* i.e. Almost Incompressible Elasticity, Stokes discretized with Q1xQ1_stabilized */
926           PetscCall(PetscObjectCompose((PetscObject)fetidp->innerbddc, "__KSPFETIDP_C", (PetscObject)C));
927         }
928         PetscCall(MatDestroy(&C));
929       }
930       /* Divergence mat */
931       if (!pcbddc->divudotp) {
932         Mat       B;
933         IS        P;
934         IS        l2l = NULL;
935         PetscBool save;
936 
937         PetscCall(PetscObjectQuery((PetscObject)fetidp->innerbddc, "__KSPFETIDP_aP", (PetscObject *)&P));
938         if (!pisz) {
939           IS       F, V;
940           PetscInt m, M;
941 
942           PetscCall(MatGetOwnershipRange(A, &m, &M));
943           PetscCall(ISCreateStride(PetscObjectComm((PetscObject)A), M - m, m, 1, &F));
944           PetscCall(ISComplement(P, m, M, &V));
945           PetscCall(MatCreateSubMatrix(A, P, V, MAT_INITIAL_MATRIX, &B));
946           {
947             Mat_IS *Bmatis = (Mat_IS *)B->data;
948             PetscCall(PetscObjectReference((PetscObject)Bmatis->getsub_cis));
949             l2l = Bmatis->getsub_cis;
950           }
951           PetscCall(ISDestroy(&V));
952           PetscCall(ISDestroy(&F));
953         } else {
954           PetscCall(MatCreateSubMatrix(A, P, NULL, MAT_INITIAL_MATRIX, &B));
955         }
956         save = pcbddc->compute_nonetflux; /* SetDivergenceMat activates nonetflux computation */
957         PetscCall(PCBDDCSetDivergenceMat(fetidp->innerbddc, B, PETSC_FALSE, l2l));
958         pcbddc->compute_nonetflux = save;
959         PetscCall(MatDestroy(&B));
960         PetscCall(ISDestroy(&l2l));
961       }
962       if (A != Ap) { /* user has provided a different Pmat, this always superseeds the setter (TODO: is it OK?) */
963         /* use monolithic operator, we restrict later */
964         PetscCall(KSPFETIDPSetPressureOperator(ksp, Ap));
965       }
966       PetscCall(PetscObjectQuery((PetscObject)fetidp->innerbddc, "__KSPFETIDP_PPmat", (PetscObject *)&PPmat));
967 
968       /* PPmat not present, use some default choice */
969       if (!PPmat) {
970         Mat C;
971 
972         PetscCall(PetscObjectQuery((PetscObject)fetidp->innerbddc, "__KSPFETIDP_C", (PetscObject *)&C));
973         if (!schp && C) { /* non-zero pressure block, most likely Almost Incompressible Elasticity */
974           PetscCall(KSPFETIDPSetPressureOperator(ksp, C));
975         } else if (!pisz && schp) { /* we need the whole pressure mass matrix to define the interface BDDC */
976           IS P;
977 
978           PetscCall(PetscObjectQuery((PetscObject)fetidp->innerbddc, "__KSPFETIDP_aP", (PetscObject *)&P));
979           PetscCall(MatCreateSubMatrix(A, P, P, MAT_INITIAL_MATRIX, &C));
980           PetscCall(MatScale(C, -1.));
981           PetscCall(KSPFETIDPSetPressureOperator(ksp, C));
982           PetscCall(MatDestroy(&C));
983         } else { /* identity (need to be scaled properly by the user using e.g. a Richardson method */
984           PetscInt nl;
985 
986           PetscCall(ISGetLocalSize(fetidp->pP, &nl));
987           PetscCall(MatCreate(PetscObjectComm((PetscObject)ksp), &C));
988           PetscCall(MatSetSizes(C, nl, nl, totP, totP));
989           PetscCall(MatSetType(C, MATAIJ));
990           PetscCall(MatMPIAIJSetPreallocation(C, 1, NULL, 0, NULL));
991           PetscCall(MatSeqAIJSetPreallocation(C, 1, NULL));
992           PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
993           PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
994           PetscCall(MatShift(C, 1.));
995           PetscCall(KSPFETIDPSetPressureOperator(ksp, C));
996           PetscCall(MatDestroy(&C));
997         }
998       }
999 
1000       /* Preconditioned operator for the pressure block */
1001       PetscCall(PetscObjectQuery((PetscObject)fetidp->innerbddc, "__KSPFETIDP_PPmat", (PetscObject *)&PPmat));
1002       if (PPmat) {
1003         Mat      C;
1004         IS       Pall;
1005         PetscInt AM, PAM, PAN, pam, pan, am, an, pl, pIl, pAg, pIg;
1006 
1007         PetscCall(PetscObjectQuery((PetscObject)fetidp->innerbddc, "__KSPFETIDP_aP", (PetscObject *)&Pall));
1008         PetscCall(MatGetSize(A, &AM, NULL));
1009         PetscCall(MatGetSize(PPmat, &PAM, &PAN));
1010         PetscCall(ISGetSize(Pall, &pAg));
1011         PetscCall(ISGetSize(fetidp->pP, &pIg));
1012         PetscCall(MatGetLocalSize(PPmat, &pam, &pan));
1013         PetscCall(MatGetLocalSize(A, &am, &an));
1014         PetscCall(ISGetLocalSize(Pall, &pIl));
1015         PetscCall(ISGetLocalSize(fetidp->pP, &pl));
1016         PetscCheck(PAM == PAN, PetscObjectComm((PetscObject)ksp), PETSC_ERR_USER, "Pressure matrix must be square, unsupported %" PetscInt_FMT " x %" PetscInt_FMT, PAM, PAN);
1017         PetscCheck(pam == pan, PetscObjectComm((PetscObject)ksp), PETSC_ERR_USER, "Local sizes of pressure matrix must be equal, unsupported %" PetscInt_FMT " x %" PetscInt_FMT, pam, pan);
1018         PetscCheck(pam == am || pam == pl || pam == pIl, PETSC_COMM_SELF, PETSC_ERR_USER, "Invalid number of local rows %" PetscInt_FMT " for pressure matrix! Supported are %" PetscInt_FMT ", %" PetscInt_FMT " or %" PetscInt_FMT, pam, am, pl, pIl);
1019         PetscCheck(pan == an || pan == pl || pan == pIl, PETSC_COMM_SELF, PETSC_ERR_USER, "Invalid number of local columns %" PetscInt_FMT " for pressure matrix! Supported are %" PetscInt_FMT ", %" PetscInt_FMT " or %" PetscInt_FMT, pan, an, pl, pIl);
1020         if (PAM == AM) { /* monolithic ordering, restrict to pressure */
1021           if (schp) {
1022             PetscCall(MatCreateSubMatrix(PPmat, Pall, Pall, MAT_INITIAL_MATRIX, &C));
1023           } else {
1024             PetscCall(MatCreateSubMatrix(PPmat, fetidp->pP, fetidp->pP, MAT_INITIAL_MATRIX, &C));
1025           }
1026         } else if (pAg == PAM) { /* global ordering for pressure only */
1027           if (!allp && !schp) {  /* solving for interface pressure only */
1028             IS restr;
1029 
1030             PetscCall(ISRenumber(fetidp->pP, NULL, NULL, &restr));
1031             PetscCall(MatCreateSubMatrix(PPmat, restr, restr, MAT_INITIAL_MATRIX, &C));
1032             PetscCall(ISDestroy(&restr));
1033           } else {
1034             PetscCall(PetscObjectReference((PetscObject)PPmat));
1035             C = PPmat;
1036           }
1037         } else if (pIg == PAM) { /* global ordering for selected pressure only */
1038           PetscCheck(!schp, PetscObjectComm((PetscObject)ksp), PETSC_ERR_PLIB, "Need the entire matrix");
1039           PetscCall(PetscObjectReference((PetscObject)PPmat));
1040           C = PPmat;
1041         } else SETERRQ(PetscObjectComm((PetscObject)ksp), PETSC_ERR_USER, "Unable to use the pressure matrix");
1042 
1043         PetscCall(KSPFETIDPSetPressureOperator(ksp, C));
1044         PetscCall(MatDestroy(&C));
1045       } else SETERRQ(PetscObjectComm((PetscObject)ksp), PETSC_ERR_PLIB, "Missing Pmat for pressure block");
1046     } else { /* totP == 0 */
1047       PetscCall(PetscObjectCompose((PetscObject)fetidp->innerbddc, "__KSPFETIDP_pP", NULL));
1048     }
1049   }
1050   PetscFunctionReturn(0);
1051 }
1052 
1053 static PetscErrorCode KSPSetUp_FETIDP(KSP ksp) {
1054   KSP_FETIDP *fetidp = (KSP_FETIDP *)ksp->data;
1055   PC_BDDC    *pcbddc = (PC_BDDC *)fetidp->innerbddc->data;
1056   PetscBool   flg;
1057 
1058   PetscFunctionBegin;
1059   PetscCall(KSPFETIDPSetUpOperators(ksp));
1060   /* set up BDDC */
1061   PetscCall(PCSetErrorIfFailure(fetidp->innerbddc, ksp->errorifnotconverged));
1062   PetscCall(PCSetUp(fetidp->innerbddc));
1063   /* FETI-DP as it is implemented needs an exact coarse solver */
1064   if (pcbddc->coarse_ksp) {
1065     PetscCall(KSPSetTolerances(pcbddc->coarse_ksp, PETSC_SMALL, PETSC_SMALL, PETSC_DEFAULT, 1000));
1066     PetscCall(KSPSetNormType(pcbddc->coarse_ksp, KSP_NORM_DEFAULT));
1067   }
1068   /* FETI-DP as it is implemented needs exact local Neumann solvers */
1069   PetscCall(KSPSetTolerances(pcbddc->ksp_R, PETSC_SMALL, PETSC_SMALL, PETSC_DEFAULT, 1000));
1070   PetscCall(KSPSetNormType(pcbddc->ksp_R, KSP_NORM_DEFAULT));
1071 
1072   /* setup FETI-DP operators
1073      If fetidp->statechanged is true, we need to update the operators
1074      needed in the saddle-point case. This should be replaced
1075      by a better logic when the FETI-DP matrix and preconditioner will
1076      have their own classes */
1077   if (pcbddc->new_primal_space || fetidp->statechanged) {
1078     Mat F; /* the FETI-DP matrix */
1079     PC  D; /* the FETI-DP preconditioner */
1080     PetscCall(KSPReset(fetidp->innerksp));
1081     PetscCall(PCBDDCCreateFETIDPOperators(fetidp->innerbddc, fetidp->fully_redundant, ((PetscObject)ksp)->prefix, &F, &D));
1082     PetscCall(KSPSetOperators(fetidp->innerksp, F, F));
1083     PetscCall(KSPSetTolerances(fetidp->innerksp, ksp->rtol, ksp->abstol, ksp->divtol, ksp->max_it));
1084     PetscCall(KSPSetPC(fetidp->innerksp, D));
1085     PetscCall(PetscObjectIncrementTabLevel((PetscObject)D, (PetscObject)fetidp->innerksp, 0));
1086     PetscCall(KSPSetFromOptions(fetidp->innerksp));
1087     PetscCall(MatCreateVecs(F, &(fetidp->innerksp)->vec_rhs, &(fetidp->innerksp)->vec_sol));
1088     PetscCall(MatDestroy(&F));
1089     PetscCall(PCDestroy(&D));
1090     if (fetidp->check) {
1091       PetscViewer viewer;
1092 
1093       if (!pcbddc->dbg_viewer) {
1094         viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)ksp));
1095       } else {
1096         viewer = pcbddc->dbg_viewer;
1097       }
1098       PetscCall(KSPFETIDPCheckOperators(ksp, viewer));
1099     }
1100   }
1101   fetidp->statechanged     = PETSC_FALSE;
1102   pcbddc->new_primal_space = PETSC_FALSE;
1103 
1104   /* propagate settings to the inner solve */
1105   PetscCall(KSPGetComputeSingularValues(ksp, &flg));
1106   PetscCall(KSPSetComputeSingularValues(fetidp->innerksp, flg));
1107   if (ksp->res_hist) PetscCall(KSPSetResidualHistory(fetidp->innerksp, ksp->res_hist, ksp->res_hist_max, ksp->res_hist_reset));
1108   PetscCall(KSPSetErrorIfNotConverged(fetidp->innerksp, ksp->errorifnotconverged));
1109   PetscCall(KSPSetUp(fetidp->innerksp));
1110   PetscFunctionReturn(0);
1111 }
1112 
1113 static PetscErrorCode KSPSolve_FETIDP(KSP ksp) {
1114   Mat                F, A;
1115   MatNullSpace       nsp;
1116   Vec                X, B, Xl, Bl;
1117   KSP_FETIDP        *fetidp = (KSP_FETIDP *)ksp->data;
1118   PC_BDDC           *pcbddc = (PC_BDDC *)fetidp->innerbddc->data;
1119   KSPConvergedReason reason;
1120   PC                 pc;
1121   PCFailedReason     pcreason;
1122   PetscInt           hist_len;
1123 
1124   PetscFunctionBegin;
1125   PetscCall(PetscCitationsRegister(citation, &cited));
1126   if (fetidp->saddlepoint) { PetscCall(PetscCitationsRegister(citation2, &cited2)); }
1127   PetscCall(KSPGetOperators(ksp, &A, NULL));
1128   PetscCall(KSPGetRhs(ksp, &B));
1129   PetscCall(KSPGetSolution(ksp, &X));
1130   PetscCall(KSPGetOperators(fetidp->innerksp, &F, NULL));
1131   PetscCall(KSPGetRhs(fetidp->innerksp, &Bl));
1132   PetscCall(KSPGetSolution(fetidp->innerksp, &Xl));
1133   PetscCall(PCBDDCMatFETIDPGetRHS(F, B, Bl));
1134   if (ksp->transpose_solve) {
1135     PetscCall(KSPSolveTranspose(fetidp->innerksp, Bl, Xl));
1136   } else {
1137     PetscCall(KSPSolve(fetidp->innerksp, Bl, Xl));
1138   }
1139   PetscCall(KSPGetConvergedReason(fetidp->innerksp, &reason));
1140   PetscCall(KSPGetPC(fetidp->innerksp, &pc));
1141   PetscCall(PCGetFailedReason(pc, &pcreason));
1142   if ((reason < 0 && reason != KSP_DIVERGED_ITS) || pcreason) {
1143     PetscInt its;
1144     PetscCall(KSPGetIterationNumber(fetidp->innerksp, &its));
1145     ksp->reason = KSP_DIVERGED_PC_FAILED;
1146     PetscCall(VecSetInf(Xl));
1147     PetscCall(PetscInfo(ksp, "Inner KSP solve failed: %s %s at iteration %" PetscInt_FMT, KSPConvergedReasons[reason], PCFailedReasons[pcreason], its));
1148   }
1149   PetscCall(PCBDDCMatFETIDPGetSolution(F, Xl, X));
1150   PetscCall(MatGetNullSpace(A, &nsp));
1151   if (nsp) PetscCall(MatNullSpaceRemove(nsp, X));
1152   /* update ksp with stats from inner ksp */
1153   PetscCall(KSPGetConvergedReason(fetidp->innerksp, &ksp->reason));
1154   PetscCall(KSPGetIterationNumber(fetidp->innerksp, &ksp->its));
1155   ksp->totalits += ksp->its;
1156   PetscCall(KSPGetResidualHistory(fetidp->innerksp, NULL, &hist_len));
1157   ksp->res_hist_len                 = (size_t)hist_len;
1158   /* restore defaults for inner BDDC (Pre/PostSolve flags) */
1159   pcbddc->temp_solution_used        = PETSC_FALSE;
1160   pcbddc->rhs_change                = PETSC_FALSE;
1161   pcbddc->exact_dirichlet_trick_app = PETSC_FALSE;
1162   PetscFunctionReturn(0);
1163 }
1164 
1165 static PetscErrorCode KSPReset_FETIDP(KSP ksp) {
1166   KSP_FETIDP *fetidp = (KSP_FETIDP *)ksp->data;
1167   PC_BDDC    *pcbddc;
1168 
1169   PetscFunctionBegin;
1170   PetscCall(ISDestroy(&fetidp->pP));
1171   PetscCall(VecDestroy(&fetidp->rhs_flip));
1172   /* avoid PCReset that does not take into account ref counting */
1173   PetscCall(PCDestroy(&fetidp->innerbddc));
1174   PetscCall(PCCreate(PetscObjectComm((PetscObject)ksp), &fetidp->innerbddc));
1175   PetscCall(PCSetType(fetidp->innerbddc, PCBDDC));
1176   pcbddc                   = (PC_BDDC *)fetidp->innerbddc->data;
1177   pcbddc->symmetric_primal = PETSC_FALSE;
1178   PetscCall(PetscLogObjectParent((PetscObject)ksp, (PetscObject)fetidp->innerbddc));
1179   PetscCall(KSPDestroy(&fetidp->innerksp));
1180   fetidp->saddlepoint  = PETSC_FALSE;
1181   fetidp->matstate     = -1;
1182   fetidp->matnnzstate  = -1;
1183   fetidp->statechanged = PETSC_TRUE;
1184   PetscFunctionReturn(0);
1185 }
1186 
1187 static PetscErrorCode KSPDestroy_FETIDP(KSP ksp) {
1188   KSP_FETIDP *fetidp = (KSP_FETIDP *)ksp->data;
1189 
1190   PetscFunctionBegin;
1191   PetscCall(KSPReset_FETIDP(ksp));
1192   PetscCall(PCDestroy(&fetidp->innerbddc));
1193   PetscCall(KSPDestroy(&fetidp->innerksp));
1194   PetscCall(PetscFree(fetidp->monctx));
1195   PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPFETIDPSetInnerBDDC_C", NULL));
1196   PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPFETIDPGetInnerBDDC_C", NULL));
1197   PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPFETIDPGetInnerKSP_C", NULL));
1198   PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPFETIDPSetPressureOperator_C", NULL));
1199   PetscCall(PetscFree(ksp->data));
1200   PetscFunctionReturn(0);
1201 }
1202 
1203 static PetscErrorCode KSPView_FETIDP(KSP ksp, PetscViewer viewer) {
1204   KSP_FETIDP *fetidp = (KSP_FETIDP *)ksp->data;
1205   PetscBool   iascii;
1206 
1207   PetscFunctionBegin;
1208   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
1209   if (iascii) {
1210     PetscCall(PetscViewerASCIIPrintf(viewer, "  fully redundant: %d\n", fetidp->fully_redundant));
1211     PetscCall(PetscViewerASCIIPrintf(viewer, "  saddle point:    %d\n", fetidp->saddlepoint));
1212     PetscCall(PetscViewerASCIIPrintf(viewer, "Inner KSP solver details\n"));
1213   }
1214   PetscCall(PetscViewerASCIIPushTab(viewer));
1215   PetscCall(KSPView(fetidp->innerksp, viewer));
1216   PetscCall(PetscViewerASCIIPopTab(viewer));
1217   if (iascii) { PetscCall(PetscViewerASCIIPrintf(viewer, "Inner BDDC solver details\n")); }
1218   PetscCall(PetscViewerASCIIPushTab(viewer));
1219   PetscCall(PCView(fetidp->innerbddc, viewer));
1220   PetscCall(PetscViewerASCIIPopTab(viewer));
1221   PetscFunctionReturn(0);
1222 }
1223 
1224 static PetscErrorCode KSPSetFromOptions_FETIDP(KSP ksp, PetscOptionItems *PetscOptionsObject) {
1225   KSP_FETIDP *fetidp = (KSP_FETIDP *)ksp->data;
1226 
1227   PetscFunctionBegin;
1228   /* set options prefixes for the inner objects, since the parent prefix will be valid at this point */
1229   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)fetidp->innerksp, ((PetscObject)ksp)->prefix));
1230   PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)fetidp->innerksp, "fetidp_"));
1231   if (!fetidp->userbddc) {
1232     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)fetidp->innerbddc, ((PetscObject)ksp)->prefix));
1233     PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)fetidp->innerbddc, "fetidp_bddc_"));
1234   }
1235   PetscOptionsHeadBegin(PetscOptionsObject, "KSP FETIDP options");
1236   PetscCall(PetscOptionsBool("-ksp_fetidp_fullyredundant", "Use fully redundant multipliers", "none", fetidp->fully_redundant, &fetidp->fully_redundant, NULL));
1237   PetscCall(PetscOptionsBool("-ksp_fetidp_saddlepoint", "Activates support for saddle-point problems", NULL, fetidp->saddlepoint, &fetidp->saddlepoint, NULL));
1238   PetscCall(PetscOptionsBool("-ksp_fetidp_check", "Activates verbose debugging output FETI-DP operators", NULL, fetidp->check, &fetidp->check, NULL));
1239   PetscOptionsHeadEnd();
1240   PetscCall(PCSetFromOptions(fetidp->innerbddc));
1241   PetscFunctionReturn(0);
1242 }
1243 
1244 /*MC
1245      KSPFETIDP - The FETI-DP method
1246 
1247    This class implements the FETI-DP method [1].
1248    The matrix for the KSP must be of type MATIS.
1249    The FETI-DP linear system (automatically generated constructing an internal PCBDDC object) is solved using an internal KSP object.
1250 
1251    Options Database Keys:
1252 +   -ksp_fetidp_fullyredundant <false>   - use a fully redundant set of Lagrange multipliers
1253 .   -ksp_fetidp_saddlepoint <false>      - activates support for saddle point problems, see [2]
1254 .   -ksp_fetidp_saddlepoint_flip <false> - usually, an incompressible Stokes problem is written as
1255                                            | A B^T | | v | = | f |
1256                                            | B 0   | | p | = | g |
1257                                            with B representing -\int_\Omega \nabla \cdot u q.
1258                                            If -ksp_fetidp_saddlepoint_flip is true, the code assumes that the user provides it as
1259                                            | A B^T | | v | = | f |
1260                                            |-B 0   | | p | = |-g |
1261 .   -ksp_fetidp_pressure_field <-1>      - activates support for saddle point problems, and identifies the pressure field id.
1262                                            If this information is not provided, the pressure field is detected by using MatFindZeroDiagonals().
1263 -   -ksp_fetidp_pressure_all <false>     - if false, uses the interface pressures, as described in [2]. If true, uses the entire pressure field.
1264 
1265    Level: Advanced
1266 
1267    Notes:
1268     Options for the inner KSP and for the customization of the PCBDDC object can be specified at command line by using the prefixes -fetidp_ and -fetidp_bddc_. E.g.,
1269 .vb
1270       -fetidp_ksp_type gmres -fetidp_bddc_pc_bddc_symmetric false
1271 .ve
1272    will use GMRES for the solution of the linear system on the Lagrange multipliers, generated using a non-symmetric PCBDDC.
1273 
1274    For saddle point problems with continuous pressures, the preconditioned operator for the pressure solver can be specified with KSPFETIDPSetPressureOperator().
1275    Alternatively, the pressure operator is extracted from the precondioned matrix (if it is different from the linear solver matrix).
1276    If none of the above, an identity matrix will be created; the user then needs to scale it through a Richardson solver.
1277    Options for the pressure solver can be prefixed with -fetidp_fielsplit_p_, E.g.
1278 .vb
1279       -fetidp_fielsplit_p_ksp_type preonly -fetidp_fielsplit_p_pc_type lu -fetidp_fielsplit_p_pc_factor_mat_solver_type mumps
1280 .ve
1281    In order to use the deluxe version of FETI-DP, you must customize the inner BDDC operator with -fetidp_bddc_pc_bddc_use_deluxe_scaling -fetidp_bddc_pc_bddc_deluxe_singlemat and use
1282    non-redundant multipliers, i.e. -ksp_fetidp_fullyredundant false. Options for the scaling solver are prefixed by -fetidp_bddelta_, E.g.
1283 .vb
1284       -fetidp_bddelta_pc_factor_mat_solver_type mumps -fetidp_bddelta_pc_type lu
1285 .ve
1286 
1287    Some of the basic options such as the maximum number of iterations and tolerances are automatically passed from this KSP to the inner KSP that actually performs the iterations.
1288 
1289    The converged reason and number of iterations computed are passed from the inner KSP to this KSP at the end of the solution.
1290 
1291    Developer Notes:
1292     Even though this method does not directly use any norms, the user is allowed to set the KSPNormType to any value.
1293     This is so users do not have to change KSPNormType options when they switch from other KSP methods to this one.
1294 
1295    References:
1296 +  * - C. Farhat, M. Lesoinne, P. LeTallec, K. Pierson, and D. Rixen, FETI-DP: a dual-primal unified FETI method. I. A faster alternative to the two-level FETI method, Internat. J. Numer. Methods Engrg., 50 (2001), pp. 1523--1544
1297 -  * - X. Tu, J. Li, A FETI-DP type domain decomposition algorithm for three-dimensional incompressible Stokes equations, SIAM J. Numer. Anal., 53 (2015), pp. 720-742
1298 
1299 .seealso: `MATIS`, `PCBDDC`, `KSPFETIDPSetInnerBDDC()`, `KSPFETIDPGetInnerBDDC()`, `KSPFETIDPGetInnerKSP()`
1300 M*/
1301 PETSC_EXTERN PetscErrorCode KSPCreate_FETIDP(KSP ksp) {
1302   KSP_FETIDP    *fetidp;
1303   KSP_FETIDPMon *monctx;
1304   PC_BDDC       *pcbddc;
1305   PC             pc;
1306 
1307   PetscFunctionBegin;
1308   PetscCall(KSPSetSupportedNorm(ksp, KSP_NORM_NONE, PC_LEFT, 3));
1309   PetscCall(KSPSetSupportedNorm(ksp, KSP_NORM_NONE, PC_RIGHT, 2));
1310   PetscCall(KSPSetSupportedNorm(ksp, KSP_NORM_PRECONDITIONED, PC_LEFT, 2));
1311   PetscCall(KSPSetSupportedNorm(ksp, KSP_NORM_PRECONDITIONED, PC_RIGHT, 2));
1312   PetscCall(KSPSetSupportedNorm(ksp, KSP_NORM_UNPRECONDITIONED, PC_LEFT, 2));
1313   PetscCall(KSPSetSupportedNorm(ksp, KSP_NORM_UNPRECONDITIONED, PC_RIGHT, 2));
1314   PetscCall(KSPSetSupportedNorm(ksp, KSP_NORM_NATURAL, PC_LEFT, 2));
1315 
1316   PetscCall(PetscNewLog(ksp, &fetidp));
1317   fetidp->matstate     = -1;
1318   fetidp->matnnzstate  = -1;
1319   fetidp->statechanged = PETSC_TRUE;
1320 
1321   ksp->data                              = (void *)fetidp;
1322   ksp->ops->setup                        = KSPSetUp_FETIDP;
1323   ksp->ops->solve                        = KSPSolve_FETIDP;
1324   ksp->ops->destroy                      = KSPDestroy_FETIDP;
1325   ksp->ops->computeeigenvalues           = KSPComputeEigenvalues_FETIDP;
1326   ksp->ops->computeextremesingularvalues = KSPComputeExtremeSingularValues_FETIDP;
1327   ksp->ops->view                         = KSPView_FETIDP;
1328   ksp->ops->setfromoptions               = KSPSetFromOptions_FETIDP;
1329   ksp->ops->buildsolution                = KSPBuildSolution_FETIDP;
1330   ksp->ops->buildresidual                = KSPBuildResidualDefault;
1331   PetscCall(KSPGetPC(ksp, &pc));
1332   PetscCall(PCSetType(pc, PCNONE));
1333   /* create the inner KSP for the Lagrange multipliers */
1334   PetscCall(KSPCreate(PetscObjectComm((PetscObject)ksp), &fetidp->innerksp));
1335   PetscCall(KSPGetPC(fetidp->innerksp, &pc));
1336   PetscCall(PCSetType(pc, PCNONE));
1337   PetscCall(PetscLogObjectParent((PetscObject)ksp, (PetscObject)fetidp->innerksp));
1338   /* monitor */
1339   PetscCall(PetscNew(&monctx));
1340   monctx->parentksp = ksp;
1341   fetidp->monctx    = monctx;
1342   PetscCall(KSPMonitorSet(fetidp->innerksp, KSPMonitor_FETIDP, fetidp->monctx, NULL));
1343   /* create the inner BDDC */
1344   PetscCall(PCCreate(PetscObjectComm((PetscObject)ksp), &fetidp->innerbddc));
1345   PetscCall(PCSetType(fetidp->innerbddc, PCBDDC));
1346   /* make sure we always obtain a consistent FETI-DP matrix
1347      for symmetric problems, the user can always customize it through the command line */
1348   pcbddc                   = (PC_BDDC *)fetidp->innerbddc->data;
1349   pcbddc->symmetric_primal = PETSC_FALSE;
1350   PetscCall(PetscLogObjectParent((PetscObject)ksp, (PetscObject)fetidp->innerbddc));
1351   /* composed functions */
1352   PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPFETIDPSetInnerBDDC_C", KSPFETIDPSetInnerBDDC_FETIDP));
1353   PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPFETIDPGetInnerBDDC_C", KSPFETIDPGetInnerBDDC_FETIDP));
1354   PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPFETIDPGetInnerKSP_C", KSPFETIDPGetInnerKSP_FETIDP));
1355   PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPFETIDPSetPressureOperator_C", KSPFETIDPSetPressureOperator_FETIDP));
1356   /* need to call KSPSetUp_FETIDP even with KSP_SETUP_NEWMATRIX */
1357   ksp->setupnewmatrix = PETSC_TRUE;
1358   PetscFunctionReturn(0);
1359 }
1360