xref: /petsc/src/ksp/pc/impls/bddc/bddcschurs.c (revision bcee047adeeb73090d7e36cc71e39fc287cdbb97)
1 #include <petsc/private/pcbddcimpl.h>
2 #include <petsc/private/pcbddcprivateimpl.h>
3 #include <../src/mat/impls/dense/seq/dense.h>
4 #include <petscblaslapack.h>
5 
6 static inline PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt *, PetscInt, PetscBT, PetscInt *, PetscInt *, PetscInt *);
7 static PetscErrorCode        PCBDDCComputeExplicitSchur(Mat, PetscBool, MatReuse, Mat *);
8 static PetscErrorCode        PCBDDCReuseSolvers_Interior(PC, Vec, Vec);
9 static PetscErrorCode        PCBDDCReuseSolvers_Correction(PC, Vec, Vec);
10 
11 /* if v2 is not present, correction is done in-place */
12 PetscErrorCode PCBDDCReuseSolversBenignAdapt(PCBDDCReuseSolvers ctx, Vec v, Vec v2, PetscBool sol, PetscBool full)
13 {
14   PetscScalar *array;
15   PetscScalar *array2;
16 
17   PetscFunctionBegin;
18   if (!ctx->benign_n) PetscFunctionReturn(PETSC_SUCCESS);
19   if (sol && full) {
20     PetscInt n_I, size_schur;
21 
22     /* get sizes */
23     PetscCall(MatGetSize(ctx->benign_csAIB, &size_schur, NULL));
24     PetscCall(VecGetSize(v, &n_I));
25     n_I = n_I - size_schur;
26     /* get schur sol from array */
27     PetscCall(VecGetArray(v, &array));
28     PetscCall(VecPlaceArray(ctx->benign_dummy_schur_vec, array + n_I));
29     PetscCall(VecRestoreArray(v, &array));
30     /* apply interior sol correction */
31     PetscCall(MatMultTranspose(ctx->benign_csAIB, ctx->benign_dummy_schur_vec, ctx->benign_corr_work));
32     PetscCall(VecResetArray(ctx->benign_dummy_schur_vec));
33     PetscCall(MatMultAdd(ctx->benign_AIIm1ones, ctx->benign_corr_work, v, v));
34   }
35   if (v2) {
36     PetscInt nl;
37 
38     PetscCall(VecGetArrayRead(v, (const PetscScalar **)&array));
39     PetscCall(VecGetLocalSize(v2, &nl));
40     PetscCall(VecGetArray(v2, &array2));
41     PetscCall(PetscArraycpy(array2, array, nl));
42   } else {
43     PetscCall(VecGetArray(v, &array));
44     array2 = array;
45   }
46   if (!sol) { /* change rhs */
47     PetscInt n;
48     for (n = 0; n < ctx->benign_n; n++) {
49       PetscScalar     sum = 0.;
50       const PetscInt *cols;
51       PetscInt        nz, i;
52 
53       PetscCall(ISGetLocalSize(ctx->benign_zerodiag_subs[n], &nz));
54       PetscCall(ISGetIndices(ctx->benign_zerodiag_subs[n], &cols));
55       for (i = 0; i < nz - 1; i++) sum += array[cols[i]];
56 #if defined(PETSC_USE_COMPLEX)
57       sum = -(PetscRealPart(sum) / nz + PETSC_i * (PetscImaginaryPart(sum) / nz));
58 #else
59       sum = -sum / nz;
60 #endif
61       for (i = 0; i < nz - 1; i++) array2[cols[i]] += sum;
62       ctx->benign_save_vals[n] = array2[cols[nz - 1]];
63       array2[cols[nz - 1]]     = sum;
64       PetscCall(ISRestoreIndices(ctx->benign_zerodiag_subs[n], &cols));
65     }
66   } else {
67     PetscInt n;
68     for (n = 0; n < ctx->benign_n; n++) {
69       PetscScalar     sum = 0.;
70       const PetscInt *cols;
71       PetscInt        nz, i;
72       PetscCall(ISGetLocalSize(ctx->benign_zerodiag_subs[n], &nz));
73       PetscCall(ISGetIndices(ctx->benign_zerodiag_subs[n], &cols));
74       for (i = 0; i < nz - 1; i++) sum += array[cols[i]];
75 #if defined(PETSC_USE_COMPLEX)
76       sum = -(PetscRealPart(sum) / nz + PETSC_i * (PetscImaginaryPart(sum) / nz));
77 #else
78       sum = -sum / nz;
79 #endif
80       for (i = 0; i < nz - 1; i++) array2[cols[i]] += sum;
81       array2[cols[nz - 1]] = ctx->benign_save_vals[n];
82       PetscCall(ISRestoreIndices(ctx->benign_zerodiag_subs[n], &cols));
83     }
84   }
85   if (v2) {
86     PetscCall(VecRestoreArrayRead(v, (const PetscScalar **)&array));
87     PetscCall(VecRestoreArray(v2, &array2));
88   } else {
89     PetscCall(VecRestoreArray(v, &array));
90   }
91   if (!sol && full) {
92     Vec      usedv;
93     PetscInt n_I, size_schur;
94 
95     /* get sizes */
96     PetscCall(MatGetSize(ctx->benign_csAIB, &size_schur, NULL));
97     PetscCall(VecGetSize(v, &n_I));
98     n_I = n_I - size_schur;
99     /* compute schur rhs correction */
100     if (v2) {
101       usedv = v2;
102     } else {
103       usedv = v;
104     }
105     /* apply schur rhs correction */
106     PetscCall(MatMultTranspose(ctx->benign_AIIm1ones, usedv, ctx->benign_corr_work));
107     PetscCall(VecGetArrayRead(usedv, (const PetscScalar **)&array));
108     PetscCall(VecPlaceArray(ctx->benign_dummy_schur_vec, array + n_I));
109     PetscCall(VecRestoreArrayRead(usedv, (const PetscScalar **)&array));
110     PetscCall(MatMultAdd(ctx->benign_csAIB, ctx->benign_corr_work, ctx->benign_dummy_schur_vec, ctx->benign_dummy_schur_vec));
111     PetscCall(VecResetArray(ctx->benign_dummy_schur_vec));
112   }
113   PetscFunctionReturn(PETSC_SUCCESS);
114 }
115 
116 static PetscErrorCode PCBDDCReuseSolvers_Solve_Private(PC pc, Vec rhs, Vec sol, PetscBool transpose, PetscBool full)
117 {
118   PCBDDCReuseSolvers ctx;
119   PetscBool          copy = PETSC_FALSE;
120 
121   PetscFunctionBegin;
122   PetscCall(PCShellGetContext(pc, &ctx));
123   if (full) {
124 #if defined(PETSC_HAVE_MUMPS)
125     PetscCall(MatMumpsSetIcntl(ctx->F, 26, -1));
126 #endif
127 #if defined(PETSC_HAVE_MKL_PARDISO)
128     PetscCall(MatMkl_PardisoSetCntl(ctx->F, 70, 0));
129 #endif
130     copy = ctx->has_vertices;
131   } else { /* interior solver */
132 #if defined(PETSC_HAVE_MUMPS)
133     PetscCall(MatMumpsSetIcntl(ctx->F, 26, 0));
134 #endif
135 #if defined(PETSC_HAVE_MKL_PARDISO)
136     PetscCall(MatMkl_PardisoSetCntl(ctx->F, 70, 1));
137 #endif
138     copy = PETSC_TRUE;
139   }
140   /* copy rhs into factored matrix workspace */
141   if (copy) {
142     PetscInt     n;
143     PetscScalar *array, *array_solver;
144 
145     PetscCall(VecGetLocalSize(rhs, &n));
146     PetscCall(VecGetArrayRead(rhs, (const PetscScalar **)&array));
147     PetscCall(VecGetArray(ctx->rhs, &array_solver));
148     PetscCall(PetscArraycpy(array_solver, array, n));
149     PetscCall(VecRestoreArray(ctx->rhs, &array_solver));
150     PetscCall(VecRestoreArrayRead(rhs, (const PetscScalar **)&array));
151 
152     PetscCall(PCBDDCReuseSolversBenignAdapt(ctx, ctx->rhs, NULL, PETSC_FALSE, full));
153     if (transpose) {
154       PetscCall(MatSolveTranspose(ctx->F, ctx->rhs, ctx->sol));
155     } else {
156       PetscCall(MatSolve(ctx->F, ctx->rhs, ctx->sol));
157     }
158     PetscCall(PCBDDCReuseSolversBenignAdapt(ctx, ctx->sol, NULL, PETSC_TRUE, full));
159 
160     /* get back data to caller worskpace */
161     PetscCall(VecGetArrayRead(ctx->sol, (const PetscScalar **)&array_solver));
162     PetscCall(VecGetArray(sol, &array));
163     PetscCall(PetscArraycpy(array, array_solver, n));
164     PetscCall(VecRestoreArray(sol, &array));
165     PetscCall(VecRestoreArrayRead(ctx->sol, (const PetscScalar **)&array_solver));
166   } else {
167     if (ctx->benign_n) {
168       PetscCall(PCBDDCReuseSolversBenignAdapt(ctx, rhs, ctx->rhs, PETSC_FALSE, full));
169       if (transpose) {
170         PetscCall(MatSolveTranspose(ctx->F, ctx->rhs, sol));
171       } else {
172         PetscCall(MatSolve(ctx->F, ctx->rhs, sol));
173       }
174       PetscCall(PCBDDCReuseSolversBenignAdapt(ctx, sol, NULL, PETSC_TRUE, full));
175     } else {
176       if (transpose) {
177         PetscCall(MatSolveTranspose(ctx->F, rhs, sol));
178       } else {
179         PetscCall(MatSolve(ctx->F, rhs, sol));
180       }
181     }
182   }
183   /* restore defaults */
184 #if defined(PETSC_HAVE_MUMPS)
185   PetscCall(MatMumpsSetIcntl(ctx->F, 26, -1));
186 #endif
187 #if defined(PETSC_HAVE_MKL_PARDISO)
188   PetscCall(MatMkl_PardisoSetCntl(ctx->F, 70, 0));
189 #endif
190   PetscFunctionReturn(PETSC_SUCCESS);
191 }
192 
193 static PetscErrorCode PCBDDCReuseSolvers_Correction(PC pc, Vec rhs, Vec sol)
194 {
195   PetscFunctionBegin;
196   PetscCall(PCBDDCReuseSolvers_Solve_Private(pc, rhs, sol, PETSC_FALSE, PETSC_TRUE));
197   PetscFunctionReturn(PETSC_SUCCESS);
198 }
199 
200 static PetscErrorCode PCBDDCReuseSolvers_CorrectionTranspose(PC pc, Vec rhs, Vec sol)
201 {
202   PetscFunctionBegin;
203   PetscCall(PCBDDCReuseSolvers_Solve_Private(pc, rhs, sol, PETSC_TRUE, PETSC_TRUE));
204   PetscFunctionReturn(PETSC_SUCCESS);
205 }
206 
207 static PetscErrorCode PCBDDCReuseSolvers_Interior(PC pc, Vec rhs, Vec sol)
208 {
209   PetscFunctionBegin;
210   PetscCall(PCBDDCReuseSolvers_Solve_Private(pc, rhs, sol, PETSC_FALSE, PETSC_FALSE));
211   PetscFunctionReturn(PETSC_SUCCESS);
212 }
213 
214 static PetscErrorCode PCBDDCReuseSolvers_InteriorTranspose(PC pc, Vec rhs, Vec sol)
215 {
216   PetscFunctionBegin;
217   PetscCall(PCBDDCReuseSolvers_Solve_Private(pc, rhs, sol, PETSC_TRUE, PETSC_FALSE));
218   PetscFunctionReturn(PETSC_SUCCESS);
219 }
220 
221 static PetscErrorCode PCBDDCReuseSolvers_View(PC pc, PetscViewer viewer)
222 {
223   PCBDDCReuseSolvers ctx;
224   PetscBool          iascii;
225 
226   PetscFunctionBegin;
227   PetscCall(PCShellGetContext(pc, &ctx));
228   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
229   if (iascii) PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_INFO));
230   PetscCall(MatView(ctx->F, viewer));
231   if (iascii) PetscCall(PetscViewerPopFormat(viewer));
232   PetscFunctionReturn(PETSC_SUCCESS);
233 }
234 
235 static PetscErrorCode PCBDDCReuseSolversReset(PCBDDCReuseSolvers reuse)
236 {
237   PetscInt i;
238 
239   PetscFunctionBegin;
240   PetscCall(MatDestroy(&reuse->F));
241   PetscCall(VecDestroy(&reuse->sol));
242   PetscCall(VecDestroy(&reuse->rhs));
243   PetscCall(PCDestroy(&reuse->interior_solver));
244   PetscCall(PCDestroy(&reuse->correction_solver));
245   PetscCall(ISDestroy(&reuse->is_R));
246   PetscCall(ISDestroy(&reuse->is_B));
247   PetscCall(VecScatterDestroy(&reuse->correction_scatter_B));
248   PetscCall(VecDestroy(&reuse->sol_B));
249   PetscCall(VecDestroy(&reuse->rhs_B));
250   for (i = 0; i < reuse->benign_n; i++) PetscCall(ISDestroy(&reuse->benign_zerodiag_subs[i]));
251   PetscCall(PetscFree(reuse->benign_zerodiag_subs));
252   PetscCall(PetscFree(reuse->benign_save_vals));
253   PetscCall(MatDestroy(&reuse->benign_csAIB));
254   PetscCall(MatDestroy(&reuse->benign_AIIm1ones));
255   PetscCall(VecDestroy(&reuse->benign_corr_work));
256   PetscCall(VecDestroy(&reuse->benign_dummy_schur_vec));
257   PetscFunctionReturn(PETSC_SUCCESS);
258 }
259 
260 static PetscErrorCode PCBDDCReuseSolvers_Destroy(PC pc)
261 {
262   PCBDDCReuseSolvers ctx;
263 
264   PetscFunctionBegin;
265   PetscCall(PCShellGetContext(pc, &ctx));
266   PetscCall(PCBDDCReuseSolversReset(ctx));
267   PetscCall(PetscFree(ctx));
268   PetscCall(PCShellSetContext(pc, NULL));
269   PetscFunctionReturn(PETSC_SUCCESS);
270 }
271 
272 static PetscErrorCode PCBDDCComputeExplicitSchur(Mat M, PetscBool issym, MatReuse reuse, Mat *S)
273 {
274   Mat         B, C, D, Bd, Cd, AinvBd;
275   KSP         ksp;
276   PC          pc;
277   PetscBool   isLU, isILU, isCHOL, Bdense, Cdense;
278   PetscReal   fill = 2.0;
279   PetscInt    n_I;
280   PetscMPIInt size;
281 
282   PetscFunctionBegin;
283   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)M), &size));
284   PetscCheck(size == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Not for parallel matrices");
285   if (reuse == MAT_REUSE_MATRIX) {
286     PetscBool Sdense;
287 
288     PetscCall(PetscObjectTypeCompare((PetscObject)*S, MATSEQDENSE, &Sdense));
289     PetscCheck(Sdense, PetscObjectComm((PetscObject)M), PETSC_ERR_SUP, "S should dense");
290   }
291   PetscCall(MatSchurComplementGetSubMatrices(M, NULL, NULL, &B, &C, &D));
292   PetscCall(MatSchurComplementGetKSP(M, &ksp));
293   PetscCall(KSPGetPC(ksp, &pc));
294   PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCLU, &isLU));
295   PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCILU, &isILU));
296   PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCCHOLESKY, &isCHOL));
297   PetscCall(PetscObjectTypeCompare((PetscObject)B, MATSEQDENSE, &Bdense));
298   PetscCall(PetscObjectTypeCompare((PetscObject)C, MATSEQDENSE, &Cdense));
299   PetscCall(MatGetSize(B, &n_I, NULL));
300   if (n_I) {
301     if (!Bdense) {
302       PetscCall(MatConvert(B, MATSEQDENSE, MAT_INITIAL_MATRIX, &Bd));
303     } else {
304       Bd = B;
305     }
306 
307     if (isLU || isILU || isCHOL) {
308       Mat fact;
309       PetscCall(KSPSetUp(ksp));
310       PetscCall(PCFactorGetMatrix(pc, &fact));
311       PetscCall(MatDuplicate(Bd, MAT_DO_NOT_COPY_VALUES, &AinvBd));
312       PetscCall(MatMatSolve(fact, Bd, AinvBd));
313     } else {
314       PetscBool ex = PETSC_TRUE;
315 
316       if (ex) {
317         Mat Ainvd;
318 
319         PetscCall(PCComputeOperator(pc, MATDENSE, &Ainvd));
320         PetscCall(MatMatMult(Ainvd, Bd, MAT_INITIAL_MATRIX, fill, &AinvBd));
321         PetscCall(MatDestroy(&Ainvd));
322       } else {
323         Vec          sol, rhs;
324         PetscScalar *arrayrhs, *arraysol;
325         PetscInt     i, nrhs, n;
326 
327         PetscCall(MatDuplicate(Bd, MAT_DO_NOT_COPY_VALUES, &AinvBd));
328         PetscCall(MatGetSize(Bd, &n, &nrhs));
329         PetscCall(MatDenseGetArray(Bd, &arrayrhs));
330         PetscCall(MatDenseGetArray(AinvBd, &arraysol));
331         PetscCall(KSPGetSolution(ksp, &sol));
332         PetscCall(KSPGetRhs(ksp, &rhs));
333         for (i = 0; i < nrhs; i++) {
334           PetscCall(VecPlaceArray(rhs, arrayrhs + i * n));
335           PetscCall(VecPlaceArray(sol, arraysol + i * n));
336           PetscCall(KSPSolve(ksp, rhs, sol));
337           PetscCall(VecResetArray(rhs));
338           PetscCall(VecResetArray(sol));
339         }
340         PetscCall(MatDenseRestoreArray(Bd, &arrayrhs));
341         PetscCall(MatDenseRestoreArray(AinvBd, &arrayrhs));
342       }
343     }
344     if (!Bdense & !issym) PetscCall(MatDestroy(&Bd));
345 
346     if (!issym) {
347       if (!Cdense) {
348         PetscCall(MatConvert(C, MATSEQDENSE, MAT_INITIAL_MATRIX, &Cd));
349       } else {
350         Cd = C;
351       }
352       PetscCall(MatMatMult(Cd, AinvBd, reuse, fill, S));
353       if (!Cdense) PetscCall(MatDestroy(&Cd));
354     } else {
355       PetscCall(MatTransposeMatMult(Bd, AinvBd, reuse, fill, S));
356       if (!Bdense) PetscCall(MatDestroy(&Bd));
357     }
358     PetscCall(MatDestroy(&AinvBd));
359   }
360 
361   if (D) {
362     Mat       Dd;
363     PetscBool Ddense;
364 
365     PetscCall(PetscObjectTypeCompare((PetscObject)D, MATSEQDENSE, &Ddense));
366     if (!Ddense) {
367       PetscCall(MatConvert(D, MATSEQDENSE, MAT_INITIAL_MATRIX, &Dd));
368     } else {
369       Dd = D;
370     }
371     if (n_I) {
372       PetscCall(MatAYPX(*S, -1.0, Dd, SAME_NONZERO_PATTERN));
373     } else {
374       if (reuse == MAT_INITIAL_MATRIX) {
375         PetscCall(MatDuplicate(Dd, MAT_COPY_VALUES, S));
376       } else {
377         PetscCall(MatCopy(Dd, *S, SAME_NONZERO_PATTERN));
378       }
379     }
380     if (!Ddense) PetscCall(MatDestroy(&Dd));
381   } else {
382     PetscCall(MatScale(*S, -1.0));
383   }
384   PetscFunctionReturn(PETSC_SUCCESS);
385 }
386 
387 PetscErrorCode PCBDDCSubSchursSetUp(PCBDDCSubSchurs sub_schurs, Mat Ain, Mat Sin, PetscBool exact_schur, PetscInt xadj[], PetscInt adjncy[], PetscInt nlayers, Vec scaling, PetscBool compute_Stilda, PetscBool reuse_solvers, PetscBool benign_trick, PetscInt benign_n, PetscInt benign_p0_lidx[], IS benign_zerodiag_subs[], Mat change, IS change_primal)
388 {
389   Mat          F, A_II, A_IB, A_BI, A_BB, AE_II;
390   Mat          S_all;
391   Vec          gstash, lstash;
392   VecScatter   sstash;
393   IS           is_I, is_I_layer;
394   IS           all_subsets, all_subsets_mult, all_subsets_n;
395   PetscScalar *stasharray, *Bwork;
396   PetscInt    *nnz, *all_local_idx_N;
397   PetscInt    *auxnum1, *auxnum2;
398   PetscInt     i, subset_size, max_subset_size;
399   PetscInt     n_B, extra, local_size, global_size;
400   PetscInt     local_stash_size;
401   PetscBLASInt B_N, B_ierr, B_lwork, *pivots;
402   MPI_Comm     comm_n;
403   PetscBool    deluxe   = PETSC_TRUE;
404   PetscBool    use_potr = PETSC_FALSE, use_sytr = PETSC_FALSE;
405   PetscViewer  matl_dbg_viewer = NULL;
406   PetscBool    flg;
407 
408   PetscFunctionBegin;
409   PetscCall(MatDestroy(&sub_schurs->A));
410   PetscCall(MatDestroy(&sub_schurs->S));
411   if (Ain) {
412     PetscCall(PetscObjectReference((PetscObject)Ain));
413     sub_schurs->A = Ain;
414   }
415 
416   PetscCall(PetscObjectReference((PetscObject)Sin));
417   sub_schurs->S = Sin;
418   if (sub_schurs->schur_explicit) sub_schurs->schur_explicit = (PetscBool)(!!sub_schurs->A);
419 
420   /* preliminary checks */
421   PetscCheck(sub_schurs->schur_explicit || !compute_Stilda, PetscObjectComm((PetscObject)sub_schurs->l2gmap), PETSC_ERR_SUP, "Adaptive selection of constraints requires MUMPS and/or MKL_PARDISO");
422 
423   if (benign_trick) sub_schurs->is_posdef = PETSC_FALSE;
424 
425   /* debug (MATLAB) */
426   if (sub_schurs->debug) {
427     PetscMPIInt size, rank;
428     PetscInt    nr, *print_schurs_ranks, print_schurs = PETSC_FALSE;
429 
430     PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)sub_schurs->l2gmap), &size));
431     PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)sub_schurs->l2gmap), &rank));
432     nr = size;
433     PetscCall(PetscMalloc1(nr, &print_schurs_ranks));
434     PetscOptionsBegin(PetscObjectComm((PetscObject)sub_schurs->l2gmap), sub_schurs->prefix, "BDDC sub_schurs options", "PC");
435     PetscCall(PetscOptionsIntArray("-sub_schurs_debug_ranks", "Ranks to debug (all if the option is not used)", NULL, print_schurs_ranks, &nr, &flg));
436     if (!flg) print_schurs = PETSC_TRUE;
437     else {
438       print_schurs = PETSC_FALSE;
439       for (i = 0; i < nr; i++)
440         if (print_schurs_ranks[i] == (PetscInt)rank) {
441           print_schurs = PETSC_TRUE;
442           break;
443         }
444     }
445     PetscOptionsEnd();
446     PetscCall(PetscFree(print_schurs_ranks));
447     if (print_schurs) {
448       char filename[256];
449 
450       PetscCall(PetscSNPrintf(filename, sizeof(filename), "sub_schurs_Schur_r%d.m", PetscGlobalRank));
451       PetscCall(PetscViewerASCIIOpen(PETSC_COMM_SELF, filename, &matl_dbg_viewer));
452       PetscCall(PetscViewerPushFormat(matl_dbg_viewer, PETSC_VIEWER_ASCII_MATLAB));
453     }
454   }
455 
456   /* restrict work on active processes */
457   if (sub_schurs->restrict_comm) {
458     PetscSubcomm subcomm;
459     PetscMPIInt  color, rank;
460 
461     color = 0;
462     if (!sub_schurs->n_subs) color = 1; /* this can happen if we are in a multilevel case or if the subdomain is disconnected */
463     PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)sub_schurs->l2gmap), &rank));
464     PetscCall(PetscSubcommCreate(PetscObjectComm((PetscObject)sub_schurs->l2gmap), &subcomm));
465     PetscCall(PetscSubcommSetNumber(subcomm, 2));
466     PetscCall(PetscSubcommSetTypeGeneral(subcomm, color, rank));
467     PetscCall(PetscCommDuplicate(PetscSubcommChild(subcomm), &comm_n, NULL));
468     PetscCall(PetscSubcommDestroy(&subcomm));
469     if (!sub_schurs->n_subs) {
470       PetscCall(PetscCommDestroy(&comm_n));
471       PetscFunctionReturn(PETSC_SUCCESS);
472     }
473   } else {
474     PetscCall(PetscCommDuplicate(PetscObjectComm((PetscObject)sub_schurs->l2gmap), &comm_n, NULL));
475   }
476 
477   /* get Schur complement matrices */
478   if (!sub_schurs->schur_explicit) {
479     Mat       tA_IB, tA_BI, tA_BB;
480     PetscBool isseqsbaij;
481     PetscCall(MatSchurComplementGetSubMatrices(sub_schurs->S, &A_II, NULL, &tA_IB, &tA_BI, &tA_BB));
482     PetscCall(PetscObjectTypeCompare((PetscObject)tA_BB, MATSEQSBAIJ, &isseqsbaij));
483     if (isseqsbaij) {
484       PetscCall(MatConvert(tA_BB, MATSEQAIJ, MAT_INITIAL_MATRIX, &A_BB));
485       PetscCall(MatConvert(tA_IB, MATSEQAIJ, MAT_INITIAL_MATRIX, &A_IB));
486       PetscCall(MatConvert(tA_BI, MATSEQAIJ, MAT_INITIAL_MATRIX, &A_BI));
487     } else {
488       PetscCall(PetscObjectReference((PetscObject)tA_BB));
489       A_BB = tA_BB;
490       PetscCall(PetscObjectReference((PetscObject)tA_IB));
491       A_IB = tA_IB;
492       PetscCall(PetscObjectReference((PetscObject)tA_BI));
493       A_BI = tA_BI;
494     }
495   } else {
496     A_II = NULL;
497     A_IB = NULL;
498     A_BI = NULL;
499     A_BB = NULL;
500   }
501   S_all = NULL;
502 
503   /* determine interior problems */
504   PetscCall(ISGetLocalSize(sub_schurs->is_I, &i));
505   if (nlayers >= 0 && i) { /* Interior problems can be different from the original one */
506     PetscBT         touched;
507     const PetscInt *idx_B;
508     PetscInt        n_I, n_B, n_local_dofs, n_prev_added, j, layer, *local_numbering;
509 
510     PetscCheck(xadj, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Cannot request layering without adjacency");
511     /* get sizes */
512     PetscCall(ISGetLocalSize(sub_schurs->is_I, &n_I));
513     PetscCall(ISGetLocalSize(sub_schurs->is_B, &n_B));
514 
515     PetscCall(PetscMalloc1(n_I + n_B, &local_numbering));
516     PetscCall(PetscBTCreate(n_I + n_B, &touched));
517     PetscCall(PetscBTMemzero(n_I + n_B, touched));
518 
519     /* all boundary dofs must be skipped when adding layers */
520     PetscCall(ISGetIndices(sub_schurs->is_B, &idx_B));
521     for (j = 0; j < n_B; j++) PetscCall(PetscBTSet(touched, idx_B[j]));
522     PetscCall(PetscArraycpy(local_numbering, idx_B, n_B));
523     PetscCall(ISRestoreIndices(sub_schurs->is_B, &idx_B));
524 
525     /* add prescribed number of layers of dofs */
526     n_local_dofs = n_B;
527     n_prev_added = n_B;
528     for (layer = 0; layer < nlayers; layer++) {
529       PetscInt n_added = 0;
530       if (n_local_dofs == n_I + n_B) break;
531       PetscCheck(n_local_dofs <= n_I + n_B, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Error querying layer %" PetscInt_FMT ". Out of bound access (%" PetscInt_FMT " > %" PetscInt_FMT ")", layer, n_local_dofs, n_I + n_B);
532       PetscCall(PCBDDCAdjGetNextLayer_Private(local_numbering + n_local_dofs, n_prev_added, touched, xadj, adjncy, &n_added));
533       n_prev_added = n_added;
534       n_local_dofs += n_added;
535       if (!n_added) break;
536     }
537     PetscCall(PetscBTDestroy(&touched));
538 
539     /* IS for I layer dofs in original numbering */
540     PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)sub_schurs->is_I), n_local_dofs - n_B, local_numbering + n_B, PETSC_COPY_VALUES, &is_I_layer));
541     PetscCall(PetscFree(local_numbering));
542     PetscCall(ISSort(is_I_layer));
543     /* IS for I layer dofs in I numbering */
544     if (!sub_schurs->schur_explicit) {
545       ISLocalToGlobalMapping ItoNmap;
546       PetscCall(ISLocalToGlobalMappingCreateIS(sub_schurs->is_I, &ItoNmap));
547       PetscCall(ISGlobalToLocalMappingApplyIS(ItoNmap, IS_GTOLM_DROP, is_I_layer, &is_I));
548       PetscCall(ISLocalToGlobalMappingDestroy(&ItoNmap));
549 
550       /* II block */
551       PetscCall(MatCreateSubMatrix(A_II, is_I, is_I, MAT_INITIAL_MATRIX, &AE_II));
552     }
553   } else {
554     PetscInt n_I;
555 
556     /* IS for I dofs in original numbering */
557     PetscCall(PetscObjectReference((PetscObject)sub_schurs->is_I));
558     is_I_layer = sub_schurs->is_I;
559 
560     /* IS for I dofs in I numbering (strided 1) */
561     if (!sub_schurs->schur_explicit) {
562       PetscCall(ISGetSize(sub_schurs->is_I, &n_I));
563       PetscCall(ISCreateStride(PetscObjectComm((PetscObject)sub_schurs->is_I), n_I, 0, 1, &is_I));
564 
565       /* II block is the same */
566       PetscCall(PetscObjectReference((PetscObject)A_II));
567       AE_II = A_II;
568     }
569   }
570 
571   /* Get info on subset sizes and sum of all subsets sizes */
572   max_subset_size = 0;
573   local_size      = 0;
574   for (i = 0; i < sub_schurs->n_subs; i++) {
575     PetscCall(ISGetLocalSize(sub_schurs->is_subs[i], &subset_size));
576     max_subset_size = PetscMax(subset_size, max_subset_size);
577     local_size += subset_size;
578   }
579 
580   /* Work arrays for local indices */
581   extra = 0;
582   PetscCall(ISGetLocalSize(sub_schurs->is_B, &n_B));
583   if (sub_schurs->schur_explicit && is_I_layer) PetscCall(ISGetLocalSize(is_I_layer, &extra));
584   PetscCall(PetscMalloc1(n_B + extra, &all_local_idx_N));
585   if (extra) {
586     const PetscInt *idxs;
587     PetscCall(ISGetIndices(is_I_layer, &idxs));
588     PetscCall(PetscArraycpy(all_local_idx_N, idxs, extra));
589     PetscCall(ISRestoreIndices(is_I_layer, &idxs));
590   }
591   PetscCall(PetscMalloc1(sub_schurs->n_subs, &auxnum1));
592   PetscCall(PetscMalloc1(sub_schurs->n_subs, &auxnum2));
593 
594   /* Get local indices in local numbering */
595   local_size       = 0;
596   local_stash_size = 0;
597   for (i = 0; i < sub_schurs->n_subs; i++) {
598     const PetscInt *idxs;
599 
600     PetscCall(ISGetLocalSize(sub_schurs->is_subs[i], &subset_size));
601     PetscCall(ISGetIndices(sub_schurs->is_subs[i], &idxs));
602     /* start (smallest in global ordering) and multiplicity */
603     auxnum1[i] = idxs[0];
604     auxnum2[i] = subset_size * subset_size;
605     /* subset indices in local numbering */
606     PetscCall(PetscArraycpy(all_local_idx_N + local_size + extra, idxs, subset_size));
607     PetscCall(ISRestoreIndices(sub_schurs->is_subs[i], &idxs));
608     local_size += subset_size;
609     local_stash_size += subset_size * subset_size;
610   }
611 
612   /* allocate extra workspace needed only for GETRI or SYTRF */
613   use_potr = use_sytr = PETSC_FALSE;
614   if (benign_trick || (sub_schurs->is_hermitian && sub_schurs->is_posdef)) {
615     use_potr = PETSC_TRUE;
616   } else if (sub_schurs->is_symmetric) {
617     use_sytr = PETSC_TRUE;
618   }
619   if (local_size && !use_potr) {
620     PetscScalar  lwork, dummyscalar = 0.;
621     PetscBLASInt dummyint = 0;
622 
623     B_lwork = -1;
624     PetscCall(PetscBLASIntCast(local_size, &B_N));
625     PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
626     if (use_sytr) {
627       PetscCallBLAS("LAPACKsytrf", LAPACKsytrf_("L", &B_N, &dummyscalar, &B_N, &dummyint, &lwork, &B_lwork, &B_ierr));
628       PetscCheck(!B_ierr, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in query to SYTRF Lapack routine %d", (int)B_ierr);
629     } else {
630       PetscCallBLAS("LAPACKgetri", LAPACKgetri_(&B_N, &dummyscalar, &B_N, &dummyint, &lwork, &B_lwork, &B_ierr));
631       PetscCheck(!B_ierr, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in query to GETRI Lapack routine %d", (int)B_ierr);
632     }
633     PetscCall(PetscFPTrapPop());
634     PetscCall(PetscBLASIntCast((PetscInt)PetscRealPart(lwork), &B_lwork));
635     PetscCall(PetscMalloc2(B_lwork, &Bwork, B_N, &pivots));
636   } else {
637     Bwork  = NULL;
638     pivots = NULL;
639   }
640 
641   /* prepare data for summing up properly schurs on subsets */
642   PetscCall(ISCreateGeneral(comm_n, sub_schurs->n_subs, auxnum1, PETSC_OWN_POINTER, &all_subsets_n));
643   PetscCall(ISLocalToGlobalMappingApplyIS(sub_schurs->l2gmap, all_subsets_n, &all_subsets));
644   PetscCall(ISDestroy(&all_subsets_n));
645   PetscCall(ISCreateGeneral(comm_n, sub_schurs->n_subs, auxnum2, PETSC_OWN_POINTER, &all_subsets_mult));
646   PetscCall(ISRenumber(all_subsets, all_subsets_mult, &global_size, &all_subsets_n));
647   PetscCall(ISDestroy(&all_subsets));
648   PetscCall(ISDestroy(&all_subsets_mult));
649   PetscCall(ISGetLocalSize(all_subsets_n, &i));
650   PetscCheck(i == local_stash_size, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid size of new subset! %" PetscInt_FMT " != %" PetscInt_FMT, i, local_stash_size);
651   PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, local_stash_size, NULL, &lstash));
652   PetscCall(VecCreateMPI(comm_n, PETSC_DECIDE, global_size, &gstash));
653   PetscCall(VecScatterCreate(lstash, NULL, gstash, all_subsets_n, &sstash));
654   PetscCall(ISDestroy(&all_subsets_n));
655 
656   /* subset indices in local boundary numbering */
657   if (!sub_schurs->is_Ej_all) {
658     PetscInt *all_local_idx_B;
659 
660     PetscCall(PetscMalloc1(local_size, &all_local_idx_B));
661     PetscCall(ISGlobalToLocalMappingApply(sub_schurs->BtoNmap, IS_GTOLM_DROP, local_size, all_local_idx_N + extra, &subset_size, all_local_idx_B));
662     PetscCheck(subset_size == local_size, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Error in sub_schurs serial (BtoNmap)! %" PetscInt_FMT " != %" PetscInt_FMT, subset_size, local_size);
663     PetscCall(ISCreateGeneral(PETSC_COMM_SELF, local_size, all_local_idx_B, PETSC_OWN_POINTER, &sub_schurs->is_Ej_all));
664   }
665 
666   if (change) {
667     ISLocalToGlobalMapping BtoS;
668     IS                     change_primal_B;
669     IS                     change_primal_all;
670 
671     PetscCheck(!sub_schurs->change_primal_sub, PETSC_COMM_SELF, PETSC_ERR_PLIB, "This should not happen");
672     PetscCheck(!sub_schurs->change, PETSC_COMM_SELF, PETSC_ERR_PLIB, "This should not happen");
673     PetscCall(PetscMalloc1(sub_schurs->n_subs, &sub_schurs->change_primal_sub));
674     for (i = 0; i < sub_schurs->n_subs; i++) {
675       ISLocalToGlobalMapping NtoS;
676       PetscCall(ISLocalToGlobalMappingCreateIS(sub_schurs->is_subs[i], &NtoS));
677       PetscCall(ISGlobalToLocalMappingApplyIS(NtoS, IS_GTOLM_DROP, change_primal, &sub_schurs->change_primal_sub[i]));
678       PetscCall(ISLocalToGlobalMappingDestroy(&NtoS));
679     }
680     PetscCall(ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap, IS_GTOLM_DROP, change_primal, &change_primal_B));
681     PetscCall(ISLocalToGlobalMappingCreateIS(sub_schurs->is_Ej_all, &BtoS));
682     PetscCall(ISGlobalToLocalMappingApplyIS(BtoS, IS_GTOLM_DROP, change_primal_B, &change_primal_all));
683     PetscCall(ISLocalToGlobalMappingDestroy(&BtoS));
684     PetscCall(ISDestroy(&change_primal_B));
685     PetscCall(PetscMalloc1(sub_schurs->n_subs, &sub_schurs->change));
686     for (i = 0; i < sub_schurs->n_subs; i++) {
687       Mat change_sub;
688 
689       PetscCall(ISGetLocalSize(sub_schurs->is_subs[i], &subset_size));
690       PetscCall(KSPCreate(PETSC_COMM_SELF, &sub_schurs->change[i]));
691       PetscCall(KSPSetType(sub_schurs->change[i], KSPPREONLY));
692       if (!sub_schurs->change_with_qr) {
693         PetscCall(MatCreateSubMatrix(change, sub_schurs->is_subs[i], sub_schurs->is_subs[i], MAT_INITIAL_MATRIX, &change_sub));
694       } else {
695         Mat change_subt;
696         PetscCall(MatCreateSubMatrix(change, sub_schurs->is_subs[i], sub_schurs->is_subs[i], MAT_INITIAL_MATRIX, &change_subt));
697         PetscCall(MatConvert(change_subt, MATSEQDENSE, MAT_INITIAL_MATRIX, &change_sub));
698         PetscCall(MatDestroy(&change_subt));
699       }
700       PetscCall(KSPSetOperators(sub_schurs->change[i], change_sub, change_sub));
701       PetscCall(MatDestroy(&change_sub));
702       PetscCall(KSPSetOptionsPrefix(sub_schurs->change[i], sub_schurs->prefix));
703       PetscCall(KSPAppendOptionsPrefix(sub_schurs->change[i], "sub_schurs_change_"));
704     }
705     PetscCall(ISDestroy(&change_primal_all));
706   }
707 
708   /* Local matrix of all local Schur on subsets (transposed) */
709   if (!sub_schurs->S_Ej_all) {
710     Mat          T;
711     PetscScalar *v;
712     PetscInt    *ii, *jj;
713     PetscInt     cum, i, j, k;
714 
715     /* MatSeqAIJSetPreallocation + MatSetValues is slow for these kind of matrices (may have large blocks)
716        Allocate properly a representative matrix and duplicate */
717     PetscCall(PetscMalloc3(local_size + 1, &ii, local_stash_size, &jj, local_stash_size, &v));
718     ii[0] = 0;
719     cum   = 0;
720     for (i = 0; i < sub_schurs->n_subs; i++) {
721       PetscCall(ISGetLocalSize(sub_schurs->is_subs[i], &subset_size));
722       for (j = 0; j < subset_size; j++) {
723         const PetscInt row = cum + j;
724         PetscInt       col = cum;
725 
726         ii[row + 1] = ii[row] + subset_size;
727         for (k = ii[row]; k < ii[row + 1]; k++) {
728           jj[k] = col;
729           col++;
730         }
731       }
732       cum += subset_size;
733     }
734     PetscCall(MatCreateSeqAIJWithArrays(PETSC_COMM_SELF, local_size, local_size, ii, jj, v, &T));
735     PetscCall(MatDuplicate(T, MAT_DO_NOT_COPY_VALUES, &sub_schurs->S_Ej_all));
736     PetscCall(MatDestroy(&T));
737     PetscCall(PetscFree3(ii, jj, v));
738   }
739   /* matrices for deluxe scaling and adaptive selection */
740   if (compute_Stilda) {
741     if (!sub_schurs->sum_S_Ej_tilda_all) PetscCall(MatDuplicate(sub_schurs->S_Ej_all, MAT_DO_NOT_COPY_VALUES, &sub_schurs->sum_S_Ej_tilda_all));
742     if (!sub_schurs->sum_S_Ej_inv_all && deluxe) PetscCall(MatDuplicate(sub_schurs->S_Ej_all, MAT_DO_NOT_COPY_VALUES, &sub_schurs->sum_S_Ej_inv_all));
743   }
744 
745   /* Compute Schur complements explicitly */
746   F = NULL;
747   if (!sub_schurs->schur_explicit) {
748     /* this code branch is used when MatFactor with Schur complement support is not present or when explicitly requested;
749        it is not efficient, unless the economic version of the scaling is used */
750     Mat          S_Ej_expl;
751     PetscScalar *work;
752     PetscInt     j, *dummy_idx;
753     PetscBool    Sdense;
754 
755     PetscCall(PetscMalloc2(max_subset_size, &dummy_idx, max_subset_size * max_subset_size, &work));
756     local_size = 0;
757     for (i = 0; i < sub_schurs->n_subs; i++) {
758       IS  is_subset_B;
759       Mat AE_EE, AE_IE, AE_EI, S_Ej;
760 
761       /* subsets in original and boundary numbering */
762       PetscCall(ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap, IS_GTOLM_DROP, sub_schurs->is_subs[i], &is_subset_B));
763       /* EE block */
764       PetscCall(MatCreateSubMatrix(A_BB, is_subset_B, is_subset_B, MAT_INITIAL_MATRIX, &AE_EE));
765       /* IE block */
766       PetscCall(MatCreateSubMatrix(A_IB, is_I, is_subset_B, MAT_INITIAL_MATRIX, &AE_IE));
767       /* EI block */
768       if (sub_schurs->is_symmetric) {
769         PetscCall(MatCreateTranspose(AE_IE, &AE_EI));
770       } else if (sub_schurs->is_hermitian) {
771         PetscCall(MatCreateHermitianTranspose(AE_IE, &AE_EI));
772       } else {
773         PetscCall(MatCreateSubMatrix(A_BI, is_subset_B, is_I, MAT_INITIAL_MATRIX, &AE_EI));
774       }
775       PetscCall(ISDestroy(&is_subset_B));
776       PetscCall(MatCreateSchurComplement(AE_II, AE_II, AE_IE, AE_EI, AE_EE, &S_Ej));
777       PetscCall(MatDestroy(&AE_EE));
778       PetscCall(MatDestroy(&AE_IE));
779       PetscCall(MatDestroy(&AE_EI));
780       if (AE_II == A_II) { /* we can reuse the same ksp */
781         KSP ksp;
782         PetscCall(MatSchurComplementGetKSP(sub_schurs->S, &ksp));
783         PetscCall(MatSchurComplementSetKSP(S_Ej, ksp));
784       } else { /* build new ksp object which inherits ksp and pc types from the original one */
785         KSP       origksp, schurksp;
786         PC        origpc, schurpc;
787         KSPType   ksp_type;
788         PetscInt  n_internal;
789         PetscBool ispcnone;
790 
791         PetscCall(MatSchurComplementGetKSP(sub_schurs->S, &origksp));
792         PetscCall(MatSchurComplementGetKSP(S_Ej, &schurksp));
793         PetscCall(KSPGetType(origksp, &ksp_type));
794         PetscCall(KSPSetType(schurksp, ksp_type));
795         PetscCall(KSPGetPC(schurksp, &schurpc));
796         PetscCall(KSPGetPC(origksp, &origpc));
797         PetscCall(PetscObjectTypeCompare((PetscObject)origpc, PCNONE, &ispcnone));
798         if (!ispcnone) {
799           PCType pc_type;
800           PetscCall(PCGetType(origpc, &pc_type));
801           PetscCall(PCSetType(schurpc, pc_type));
802         } else {
803           PetscCall(PCSetType(schurpc, PCLU));
804         }
805         PetscCall(ISGetSize(is_I, &n_internal));
806         if (!n_internal) { /* UMFPACK gives error with 0 sized problems */
807           MatSolverType solver = NULL;
808           PetscCall(PCFactorGetMatSolverType(origpc, (MatSolverType *)&solver));
809           if (solver) PetscCall(PCFactorSetMatSolverType(schurpc, solver));
810         }
811         PetscCall(KSPSetUp(schurksp));
812       }
813       PetscCall(ISGetLocalSize(sub_schurs->is_subs[i], &subset_size));
814       PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, subset_size, subset_size, work, &S_Ej_expl));
815       PetscCall(PCBDDCComputeExplicitSchur(S_Ej, sub_schurs->is_symmetric, MAT_REUSE_MATRIX, &S_Ej_expl));
816       PetscCall(PetscObjectTypeCompare((PetscObject)S_Ej_expl, MATSEQDENSE, &Sdense));
817       PetscCheck(Sdense, PETSC_COMM_SELF, PETSC_ERR_SUP, "Not yet implemented for sparse matrices");
818       for (j = 0; j < subset_size; j++) dummy_idx[j] = local_size + j;
819       PetscCall(MatSetValues(sub_schurs->S_Ej_all, subset_size, dummy_idx, subset_size, dummy_idx, work, INSERT_VALUES));
820       PetscCall(MatDestroy(&S_Ej));
821       PetscCall(MatDestroy(&S_Ej_expl));
822       local_size += subset_size;
823     }
824     PetscCall(PetscFree2(dummy_idx, work));
825     /* free */
826     PetscCall(ISDestroy(&is_I));
827     PetscCall(MatDestroy(&AE_II));
828     PetscCall(PetscFree(all_local_idx_N));
829   } else {
830     Mat                A, cs_AIB_mat = NULL, benign_AIIm1_ones_mat = NULL;
831     Mat               *gdswA;
832     Vec                Dall = NULL;
833     IS                 is_A_all, *is_p_r = NULL;
834     MatType            Stype;
835     PetscScalar       *work, *S_data, *schur_factor, infty = PETSC_MAX_REAL;
836     PetscScalar       *SEj_arr = NULL, *SEjinv_arr = NULL;
837     const PetscScalar *rS_data;
838     PetscInt           n, n_I, size_schur, size_active_schur, cum, cum2;
839     PetscBool          economic, solver_S, S_lower_triangular = PETSC_FALSE;
840     PetscBool          schur_has_vertices, factor_workaround;
841     PetscBool          use_cholesky;
842 #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
843     PetscBool oldpin;
844 #endif
845 
846     /* get sizes */
847     n_I = 0;
848     if (is_I_layer) PetscCall(ISGetLocalSize(is_I_layer, &n_I));
849     economic = PETSC_FALSE;
850     PetscCall(ISGetLocalSize(sub_schurs->is_I, &cum));
851     if (cum != n_I) economic = PETSC_TRUE;
852     PetscCall(MatGetLocalSize(sub_schurs->A, &n, NULL));
853     size_active_schur = local_size;
854 
855     /* import scaling vector (wrong formulation if we have 3D edges) */
856     if (scaling && compute_Stilda) {
857       const PetscScalar *array;
858       PetscScalar       *array2;
859       const PetscInt    *idxs;
860       PetscInt           i;
861 
862       PetscCall(ISGetIndices(sub_schurs->is_Ej_all, &idxs));
863       PetscCall(VecCreateSeq(PETSC_COMM_SELF, size_active_schur, &Dall));
864       PetscCall(VecGetArrayRead(scaling, &array));
865       PetscCall(VecGetArray(Dall, &array2));
866       for (i = 0; i < size_active_schur; i++) array2[i] = array[idxs[i]];
867       PetscCall(VecRestoreArray(Dall, &array2));
868       PetscCall(VecRestoreArrayRead(scaling, &array));
869       PetscCall(ISRestoreIndices(sub_schurs->is_Ej_all, &idxs));
870       deluxe = PETSC_FALSE;
871     }
872 
873     /* size active schurs does not count any dirichlet or vertex dof on the interface */
874     factor_workaround  = PETSC_FALSE;
875     schur_has_vertices = PETSC_FALSE;
876     cum                = n_I + size_active_schur;
877     if (sub_schurs->is_dir) {
878       const PetscInt *idxs;
879       PetscInt        n_dir;
880 
881       PetscCall(ISGetLocalSize(sub_schurs->is_dir, &n_dir));
882       PetscCall(ISGetIndices(sub_schurs->is_dir, &idxs));
883       PetscCall(PetscArraycpy(all_local_idx_N + cum, idxs, n_dir));
884       PetscCall(ISRestoreIndices(sub_schurs->is_dir, &idxs));
885       cum += n_dir;
886       if (!sub_schurs->gdsw) factor_workaround = PETSC_TRUE;
887     }
888     /* include the primal vertices in the Schur complement */
889     if (exact_schur && sub_schurs->is_vertices && (compute_Stilda || benign_n)) {
890       PetscInt n_v;
891 
892       PetscCall(ISGetLocalSize(sub_schurs->is_vertices, &n_v));
893       if (n_v) {
894         const PetscInt *idxs;
895 
896         PetscCall(ISGetIndices(sub_schurs->is_vertices, &idxs));
897         PetscCall(PetscArraycpy(all_local_idx_N + cum, idxs, n_v));
898         PetscCall(ISRestoreIndices(sub_schurs->is_vertices, &idxs));
899         cum += n_v;
900         if (!sub_schurs->gdsw) factor_workaround = PETSC_TRUE;
901         schur_has_vertices = PETSC_TRUE;
902       }
903     }
904     size_schur = cum - n_I;
905     PetscCall(ISCreateGeneral(PETSC_COMM_SELF, cum, all_local_idx_N, PETSC_OWN_POINTER, &is_A_all));
906 #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
907     oldpin = sub_schurs->A->boundtocpu;
908     PetscCall(MatBindToCPU(sub_schurs->A, PETSC_TRUE));
909 #endif
910     if (cum == n) {
911       PetscCall(ISSetPermutation(is_A_all));
912       PetscCall(MatPermute(sub_schurs->A, is_A_all, is_A_all, &A));
913     } else {
914       PetscCall(MatCreateSubMatrix(sub_schurs->A, is_A_all, is_A_all, MAT_INITIAL_MATRIX, &A));
915     }
916 #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
917     PetscCall(MatBindToCPU(sub_schurs->A, oldpin));
918 #endif
919     PetscCall(MatSetOptionsPrefixFactor(A, sub_schurs->prefix));
920     PetscCall(MatAppendOptionsPrefixFactor(A, "sub_schurs_"));
921 
922     /* if we actually change the basis for the pressures, LDL^T factors will use a lot of memory
923        this is a workaround */
924     if (benign_n) {
925       Vec                    v, benign_AIIm1_ones;
926       ISLocalToGlobalMapping N_to_reor;
927       IS                     is_p0, is_p0_p;
928       PetscScalar           *cs_AIB, *AIIm1_data;
929       PetscInt               sizeA;
930 
931       PetscCall(ISLocalToGlobalMappingCreateIS(is_A_all, &N_to_reor));
932       PetscCall(ISCreateGeneral(PETSC_COMM_SELF, benign_n, benign_p0_lidx, PETSC_COPY_VALUES, &is_p0));
933       PetscCall(ISGlobalToLocalMappingApplyIS(N_to_reor, IS_GTOLM_DROP, is_p0, &is_p0_p));
934       PetscCall(ISDestroy(&is_p0));
935       PetscCall(MatCreateVecs(A, &v, &benign_AIIm1_ones));
936       PetscCall(VecGetSize(v, &sizeA));
937       PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, sizeA, benign_n, NULL, &benign_AIIm1_ones_mat));
938       PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, size_schur, benign_n, NULL, &cs_AIB_mat));
939       PetscCall(MatDenseGetArray(cs_AIB_mat, &cs_AIB));
940       PetscCall(MatDenseGetArray(benign_AIIm1_ones_mat, &AIIm1_data));
941       PetscCall(PetscMalloc1(benign_n, &is_p_r));
942       /* compute colsum of A_IB restricted to pressures */
943       for (i = 0; i < benign_n; i++) {
944         const PetscScalar *array;
945         const PetscInt    *idxs;
946         PetscInt           j, nz;
947 
948         PetscCall(ISGlobalToLocalMappingApplyIS(N_to_reor, IS_GTOLM_DROP, benign_zerodiag_subs[i], &is_p_r[i]));
949         PetscCall(ISGetLocalSize(is_p_r[i], &nz));
950         PetscCall(ISGetIndices(is_p_r[i], &idxs));
951         for (j = 0; j < nz; j++) AIIm1_data[idxs[j] + sizeA * i] = 1.;
952         PetscCall(ISRestoreIndices(is_p_r[i], &idxs));
953         PetscCall(VecPlaceArray(benign_AIIm1_ones, AIIm1_data + sizeA * i));
954         PetscCall(MatMult(A, benign_AIIm1_ones, v));
955         PetscCall(VecResetArray(benign_AIIm1_ones));
956         PetscCall(VecGetArrayRead(v, &array));
957         for (j = 0; j < size_schur; j++) {
958 #if defined(PETSC_USE_COMPLEX)
959           cs_AIB[i * size_schur + j] = (PetscRealPart(array[j + n_I]) / nz + PETSC_i * (PetscImaginaryPart(array[j + n_I]) / nz));
960 #else
961           cs_AIB[i * size_schur + j] = array[j + n_I] / nz;
962 #endif
963         }
964         PetscCall(VecRestoreArrayRead(v, &array));
965       }
966       PetscCall(MatDenseRestoreArray(cs_AIB_mat, &cs_AIB));
967       PetscCall(MatDenseRestoreArray(benign_AIIm1_ones_mat, &AIIm1_data));
968       PetscCall(VecDestroy(&v));
969       PetscCall(VecDestroy(&benign_AIIm1_ones));
970       PetscCall(MatSetOption(A, MAT_KEEP_NONZERO_PATTERN, PETSC_FALSE));
971       PetscCall(MatSetOption(A, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_FALSE));
972       PetscCall(MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE));
973       PetscCall(MatZeroRowsColumnsIS(A, is_p0_p, 1.0, NULL, NULL));
974       PetscCall(ISDestroy(&is_p0_p));
975       PetscCall(ISLocalToGlobalMappingDestroy(&N_to_reor));
976     }
977     PetscCall(MatSetOption(A, MAT_SYMMETRIC, sub_schurs->is_symmetric));
978     PetscCall(MatSetOption(A, MAT_HERMITIAN, sub_schurs->is_hermitian));
979     PetscCall(MatSetOption(A, MAT_SPD, sub_schurs->is_posdef));
980 
981     /* for complexes, symmetric and hermitian at the same time implies null imaginary part */
982     use_cholesky = (PetscBool)((use_potr || use_sytr) && sub_schurs->is_hermitian && sub_schurs->is_symmetric);
983 
984     /* when using the benign subspace trick, the local Schur complements are SPD */
985     /* MKL_PARDISO does not handle well the computation of a Schur complement from a symmetric indefinite factorization
986        Use LU and adapt pivoting perturbation (still, solution is not as accurate as with using MUMPS) */
987     if (benign_trick) {
988       sub_schurs->is_posdef = PETSC_TRUE;
989       PetscCall(PetscStrcmp(sub_schurs->mat_solver_type, MATSOLVERMKL_PARDISO, &flg));
990       if (flg) use_cholesky = PETSC_FALSE;
991     }
992 
993     if (n_I) {
994       IS        is_schur;
995       char      stype[64];
996       PetscBool gpu = PETSC_FALSE;
997 
998       if (use_cholesky) {
999         PetscCall(MatGetFactor(A, sub_schurs->mat_solver_type, MAT_FACTOR_CHOLESKY, &F));
1000       } else {
1001         PetscCall(MatGetFactor(A, sub_schurs->mat_solver_type, MAT_FACTOR_LU, &F));
1002       }
1003       PetscCheck(F, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "MatGetFactor not supported by matrix instance of type %s. Rerun with \"-info :mat | grep MatGetFactor_\" for additional information", ((PetscObject)A)->type_name);
1004       PetscCall(MatSetErrorIfFailure(A, PETSC_TRUE));
1005 #if defined(PETSC_HAVE_MKL_PARDISO)
1006       if (benign_trick) PetscCall(MatMkl_PardisoSetCntl(F, 10, 10));
1007 #endif
1008       /* subsets ordered last */
1009       PetscCall(ISCreateStride(PETSC_COMM_SELF, size_schur, n_I, 1, &is_schur));
1010       PetscCall(MatFactorSetSchurIS(F, is_schur));
1011       PetscCall(ISDestroy(&is_schur));
1012 
1013       /* factorization step */
1014       if (use_cholesky) {
1015         PetscCall(MatCholeskyFactorSymbolic(F, A, NULL, NULL));
1016 #if defined(PETSC_HAVE_MUMPS) /* be sure that icntl 19 is not set by command line */
1017         PetscCall(MatMumpsSetIcntl(F, 19, 2));
1018 #endif
1019         PetscCall(MatCholeskyFactorNumeric(F, A, NULL));
1020         S_lower_triangular = PETSC_TRUE;
1021       } else {
1022         PetscCall(MatLUFactorSymbolic(F, A, NULL, NULL, NULL));
1023 #if defined(PETSC_HAVE_MUMPS) /* be sure that icntl 19 is not set by command line */
1024         PetscCall(MatMumpsSetIcntl(F, 19, 3));
1025 #endif
1026         PetscCall(MatLUFactorNumeric(F, A, NULL));
1027       }
1028       PetscCall(MatViewFromOptions(F, (PetscObject)A, "-mat_factor_view"));
1029 
1030       if (matl_dbg_viewer) {
1031         Mat S;
1032         IS  is;
1033 
1034         PetscCall(PetscObjectSetName((PetscObject)A, "A"));
1035         PetscCall(MatView(A, matl_dbg_viewer));
1036         PetscCall(MatFactorCreateSchurComplement(F, &S, NULL));
1037         PetscCall(PetscObjectSetName((PetscObject)S, "S"));
1038         PetscCall(MatView(S, matl_dbg_viewer));
1039         PetscCall(MatDestroy(&S));
1040         PetscCall(ISCreateStride(PETSC_COMM_SELF, n_I, 0, 1, &is));
1041         PetscCall(PetscObjectSetName((PetscObject)is, "I"));
1042         PetscCall(ISView(is, matl_dbg_viewer));
1043         PetscCall(ISDestroy(&is));
1044         PetscCall(ISCreateStride(PETSC_COMM_SELF, size_schur, n_I, 1, &is));
1045         PetscCall(PetscObjectSetName((PetscObject)is, "B"));
1046         PetscCall(ISView(is, matl_dbg_viewer));
1047         PetscCall(ISDestroy(&is));
1048         PetscCall(PetscObjectSetName((PetscObject)is_A_all, "IA"));
1049         PetscCall(ISView(is_A_all, matl_dbg_viewer));
1050         for (i = 0, cum = 0; i < sub_schurs->n_subs; i++) {
1051           IS   is;
1052           char name[16];
1053 
1054           PetscCall(PetscSNPrintf(name, sizeof(name), "IE%" PetscInt_FMT, i));
1055           PetscCall(ISGetLocalSize(sub_schurs->is_subs[i], &subset_size));
1056           PetscCall(ISCreateStride(PETSC_COMM_SELF, subset_size, cum, 1, &is));
1057           PetscCall(PetscObjectSetName((PetscObject)is, name));
1058           PetscCall(ISView(is, matl_dbg_viewer));
1059           PetscCall(ISDestroy(&is));
1060           if (sub_schurs->change) {
1061             Mat T;
1062 
1063             PetscCall(PetscSNPrintf(name, sizeof(name), "TE%" PetscInt_FMT, i));
1064             PetscCall(KSPGetOperators(sub_schurs->change[i], &T, NULL));
1065             PetscCall(PetscObjectSetName((PetscObject)T, name));
1066             PetscCall(MatView(T, matl_dbg_viewer));
1067             PetscCall(PetscSNPrintf(name, sizeof(name), "ITE%" PetscInt_FMT, i));
1068             PetscCall(PetscObjectSetName((PetscObject)sub_schurs->change_primal_sub[i], name));
1069             PetscCall(ISView(sub_schurs->change_primal_sub[i], matl_dbg_viewer));
1070           }
1071           cum += subset_size;
1072         }
1073         PetscCall(PetscViewerFlush(matl_dbg_viewer));
1074       }
1075 
1076       /* get explicit Schur Complement computed during numeric factorization */
1077       PetscCall(MatFactorGetSchurComplement(F, &S_all, NULL));
1078       PetscCall(PetscStrncpy(stype, MATSEQDENSE, sizeof(stype)));
1079 #if defined(PETSC_HAVE_CUDA)
1080       PetscCall(PetscObjectTypeCompareAny((PetscObject)A, &gpu, MATSEQAIJVIENNACL, MATSEQAIJCUSPARSE, ""));
1081 #endif
1082       if (gpu) PetscCall(PetscStrncpy(stype, MATSEQDENSECUDA, sizeof(stype)));
1083       PetscCall(PetscOptionsGetString(NULL, sub_schurs->prefix, "-sub_schurs_schur_mat_type", stype, sizeof(stype), NULL));
1084       PetscCall(MatConvert(S_all, stype, MAT_INPLACE_MATRIX, &S_all));
1085       PetscCall(MatSetOption(S_all, MAT_SPD, sub_schurs->is_posdef));
1086       PetscCall(MatSetOption(S_all, MAT_HERMITIAN, sub_schurs->is_hermitian));
1087       PetscCall(MatGetType(S_all, &Stype));
1088 
1089       /* we can reuse the solvers if we are not using the economic version */
1090       reuse_solvers = (PetscBool)(reuse_solvers && !economic);
1091       if (!sub_schurs->gdsw) {
1092         factor_workaround = (PetscBool)(reuse_solvers && factor_workaround);
1093         if (!sub_schurs->is_posdef && factor_workaround && compute_Stilda && size_active_schur) reuse_solvers = factor_workaround = PETSC_FALSE;
1094       }
1095       solver_S = PETSC_TRUE;
1096 
1097       /* update the Schur complement with the change of basis on the pressures */
1098       if (benign_n) {
1099         const PetscScalar *cs_AIB;
1100         PetscScalar       *S_data, *AIIm1_data;
1101         Mat                S2 = NULL, S3 = NULL; /* dbg */
1102         PetscScalar       *S2_data, *S3_data;    /* dbg */
1103         Vec                v, benign_AIIm1_ones;
1104         PetscInt           sizeA;
1105 
1106         PetscCall(MatDenseGetArray(S_all, &S_data));
1107         PetscCall(MatCreateVecs(A, &v, &benign_AIIm1_ones));
1108         PetscCall(VecGetSize(v, &sizeA));
1109 #if defined(PETSC_HAVE_MUMPS)
1110         PetscCall(MatMumpsSetIcntl(F, 26, 0));
1111 #endif
1112 #if defined(PETSC_HAVE_MKL_PARDISO)
1113         PetscCall(MatMkl_PardisoSetCntl(F, 70, 1));
1114 #endif
1115         PetscCall(MatDenseGetArrayRead(cs_AIB_mat, &cs_AIB));
1116         PetscCall(MatDenseGetArray(benign_AIIm1_ones_mat, &AIIm1_data));
1117         if (matl_dbg_viewer) {
1118           PetscCall(MatDuplicate(S_all, MAT_DO_NOT_COPY_VALUES, &S2));
1119           PetscCall(MatDuplicate(S_all, MAT_DO_NOT_COPY_VALUES, &S3));
1120           PetscCall(MatDenseGetArray(S2, &S2_data));
1121           PetscCall(MatDenseGetArray(S3, &S3_data));
1122         }
1123         for (i = 0; i < benign_n; i++) {
1124           PetscScalar    *array, sum = 0., one = 1., *sums;
1125           const PetscInt *idxs;
1126           PetscInt        k, j, nz;
1127           PetscBLASInt    B_k, B_n;
1128 
1129           PetscCall(PetscCalloc1(benign_n, &sums));
1130           PetscCall(VecPlaceArray(benign_AIIm1_ones, AIIm1_data + sizeA * i));
1131           PetscCall(VecCopy(benign_AIIm1_ones, v));
1132           PetscCall(MatSolve(F, v, benign_AIIm1_ones));
1133           PetscCall(MatMult(A, benign_AIIm1_ones, v));
1134           PetscCall(VecResetArray(benign_AIIm1_ones));
1135           /* p0 dofs (eliminated) are excluded from the sums */
1136           for (k = 0; k < benign_n; k++) {
1137             PetscCall(ISGetLocalSize(is_p_r[k], &nz));
1138             PetscCall(ISGetIndices(is_p_r[k], &idxs));
1139             for (j = 0; j < nz - 1; j++) sums[k] -= AIIm1_data[idxs[j] + sizeA * i];
1140             PetscCall(ISRestoreIndices(is_p_r[k], &idxs));
1141           }
1142           PetscCall(VecGetArrayRead(v, (const PetscScalar **)&array));
1143           if (matl_dbg_viewer) {
1144             Vec  vv;
1145             char name[16];
1146 
1147             PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, size_schur, array + n_I, &vv));
1148             PetscCall(PetscSNPrintf(name, sizeof(name), "Pvs%" PetscInt_FMT, i));
1149             PetscCall(PetscObjectSetName((PetscObject)vv, name));
1150             PetscCall(VecView(vv, matl_dbg_viewer));
1151           }
1152           /* perform sparse rank updates on symmetric Schur (TODO: move outside of the loop?) */
1153           /* cs_AIB already scaled by 1./nz */
1154           B_k = 1;
1155           for (k = 0; k < benign_n; k++) {
1156             sum = sums[k];
1157             PetscCall(PetscBLASIntCast(size_schur, &B_n));
1158 
1159             if (PetscAbsScalar(sum) == 0.0) continue;
1160             if (k == i) {
1161               PetscCallBLAS("BLASsyrk", BLASsyrk_("L", "N", &B_n, &B_k, &sum, cs_AIB + i * size_schur, &B_n, &one, S_data, &B_n));
1162               if (matl_dbg_viewer) PetscCallBLAS("BLASsyrk", BLASsyrk_("L", "N", &B_n, &B_k, &sum, cs_AIB + i * size_schur, &B_n, &one, S3_data, &B_n));
1163             } else { /* XXX Is it correct to use symmetric rank-2 update with half of the sum? */
1164               sum /= 2.0;
1165               PetscCallBLAS("BLASsyr2k", BLASsyr2k_("L", "N", &B_n, &B_k, &sum, cs_AIB + k * size_schur, &B_n, cs_AIB + i * size_schur, &B_n, &one, S_data, &B_n));
1166               if (matl_dbg_viewer) PetscCallBLAS("BLASsyr2k", BLASsyr2k_("L", "N", &B_n, &B_k, &sum, cs_AIB + k * size_schur, &B_n, cs_AIB + i * size_schur, &B_n, &one, S3_data, &B_n));
1167             }
1168           }
1169           sum = 1.;
1170           PetscCallBLAS("BLASsyr2k", BLASsyr2k_("L", "N", &B_n, &B_k, &sum, array + n_I, &B_n, cs_AIB + i * size_schur, &B_n, &one, S_data, &B_n));
1171           if (matl_dbg_viewer) PetscCallBLAS("BLASsyr2k", BLASsyr2k_("L", "N", &B_n, &B_k, &sum, array + n_I, &B_n, cs_AIB + i * size_schur, &B_n, &one, S2_data, &B_n));
1172           PetscCall(VecRestoreArrayRead(v, (const PetscScalar **)&array));
1173           /* set p0 entry of AIIm1_ones to zero */
1174           PetscCall(ISGetLocalSize(is_p_r[i], &nz));
1175           PetscCall(ISGetIndices(is_p_r[i], &idxs));
1176           for (j = 0; j < benign_n; j++) AIIm1_data[idxs[nz - 1] + sizeA * j] = 0.;
1177           PetscCall(ISRestoreIndices(is_p_r[i], &idxs));
1178           PetscCall(PetscFree(sums));
1179         }
1180         PetscCall(VecDestroy(&benign_AIIm1_ones));
1181         if (matl_dbg_viewer) {
1182           PetscCall(MatDenseRestoreArray(S2, &S2_data));
1183           PetscCall(MatDenseRestoreArray(S3, &S3_data));
1184         }
1185         if (!S_lower_triangular) { /* I need to expand the upper triangular data (column oriented) */
1186           PetscInt k, j;
1187           for (k = 0; k < size_schur; k++) {
1188             for (j = k; j < size_schur; j++) S_data[j * size_schur + k] = PetscConj(S_data[k * size_schur + j]);
1189           }
1190         }
1191 
1192         /* restore defaults */
1193 #if defined(PETSC_HAVE_MUMPS)
1194         PetscCall(MatMumpsSetIcntl(F, 26, -1));
1195 #endif
1196 #if defined(PETSC_HAVE_MKL_PARDISO)
1197         PetscCall(MatMkl_PardisoSetCntl(F, 70, 0));
1198 #endif
1199         PetscCall(MatDenseRestoreArrayRead(cs_AIB_mat, &cs_AIB));
1200         PetscCall(MatDenseRestoreArray(benign_AIIm1_ones_mat, &AIIm1_data));
1201         PetscCall(VecDestroy(&v));
1202         PetscCall(MatDenseRestoreArray(S_all, &S_data));
1203         if (matl_dbg_viewer) {
1204           Mat S;
1205 
1206           PetscCall(MatFactorRestoreSchurComplement(F, &S_all, MAT_FACTOR_SCHUR_UNFACTORED));
1207           PetscCall(MatFactorCreateSchurComplement(F, &S, NULL));
1208           PetscCall(PetscObjectSetName((PetscObject)S, "Sb"));
1209           PetscCall(MatView(S, matl_dbg_viewer));
1210           PetscCall(MatDestroy(&S));
1211           PetscCall(PetscObjectSetName((PetscObject)S2, "S2P"));
1212           PetscCall(MatView(S2, matl_dbg_viewer));
1213           PetscCall(PetscObjectSetName((PetscObject)S3, "S3P"));
1214           PetscCall(MatView(S3, matl_dbg_viewer));
1215           PetscCall(PetscObjectSetName((PetscObject)cs_AIB_mat, "cs"));
1216           PetscCall(MatView(cs_AIB_mat, matl_dbg_viewer));
1217           PetscCall(MatFactorGetSchurComplement(F, &S_all, NULL));
1218         }
1219         PetscCall(MatDestroy(&S2));
1220         PetscCall(MatDestroy(&S3));
1221       }
1222       if (!reuse_solvers) {
1223         for (i = 0; i < benign_n; i++) PetscCall(ISDestroy(&is_p_r[i]));
1224         PetscCall(PetscFree(is_p_r));
1225         PetscCall(MatDestroy(&cs_AIB_mat));
1226         PetscCall(MatDestroy(&benign_AIIm1_ones_mat));
1227       }
1228     } else { /* we can't use MatFactor when size_schur == size_of_the_problem */
1229       PetscCall(MatConvert(A, MATSEQDENSE, MAT_INITIAL_MATRIX, &S_all));
1230       PetscCall(MatGetType(S_all, &Stype));
1231       reuse_solvers     = PETSC_FALSE; /* TODO: why we can't reuse the solvers here? */
1232       factor_workaround = PETSC_FALSE;
1233       solver_S          = PETSC_FALSE;
1234     }
1235 
1236     if (reuse_solvers) {
1237       Mat                A_II, pA_II, Afake;
1238       Vec                vec1_B;
1239       PCBDDCReuseSolvers msolv_ctx;
1240       PetscInt           n_R;
1241 
1242       if (sub_schurs->reuse_solver) {
1243         PetscCall(PCBDDCReuseSolversReset(sub_schurs->reuse_solver));
1244       } else {
1245         PetscCall(PetscNew(&sub_schurs->reuse_solver));
1246       }
1247       msolv_ctx = sub_schurs->reuse_solver;
1248       PetscCall(MatSchurComplementGetSubMatrices(sub_schurs->S, &A_II, &pA_II, NULL, NULL, NULL));
1249       PetscCall(PetscObjectReference((PetscObject)F));
1250       msolv_ctx->F = F;
1251       PetscCall(MatCreateVecs(F, &msolv_ctx->sol, NULL));
1252       /* currently PETSc has no support for MatSolve(F,x,x), so cheat and let rhs and sol share the same memory */
1253       {
1254         PetscScalar *array;
1255         PetscInt     n;
1256 
1257         PetscCall(VecGetLocalSize(msolv_ctx->sol, &n));
1258         PetscCall(VecGetArray(msolv_ctx->sol, &array));
1259         PetscCall(VecCreateSeqWithArray(PetscObjectComm((PetscObject)msolv_ctx->sol), 1, n, array, &msolv_ctx->rhs));
1260         PetscCall(VecRestoreArray(msolv_ctx->sol, &array));
1261       }
1262       msolv_ctx->has_vertices = schur_has_vertices;
1263 
1264       /* interior solver */
1265       PetscCall(PCCreate(PetscObjectComm((PetscObject)A_II), &msolv_ctx->interior_solver));
1266       PetscCall(PCSetOperators(msolv_ctx->interior_solver, A_II, pA_II));
1267       PetscCall(PCSetType(msolv_ctx->interior_solver, PCSHELL));
1268       PetscCall(PCShellSetName(msolv_ctx->interior_solver, "Interior solver (w/o Schur factorization)"));
1269       PetscCall(PCShellSetContext(msolv_ctx->interior_solver, msolv_ctx));
1270       PetscCall(PCShellSetView(msolv_ctx->interior_solver, PCBDDCReuseSolvers_View));
1271       PetscCall(PCShellSetApply(msolv_ctx->interior_solver, PCBDDCReuseSolvers_Interior));
1272       PetscCall(PCShellSetApplyTranspose(msolv_ctx->interior_solver, PCBDDCReuseSolvers_InteriorTranspose));
1273       if (sub_schurs->gdsw) PetscCall(PCShellSetDestroy(msolv_ctx->interior_solver, PCBDDCReuseSolvers_Destroy));
1274 
1275       /* correction solver */
1276       if (!sub_schurs->gdsw) {
1277         PetscCall(PCCreate(PetscObjectComm((PetscObject)A_II), &msolv_ctx->correction_solver));
1278         PetscCall(PCSetType(msolv_ctx->correction_solver, PCSHELL));
1279         PetscCall(PCShellSetName(msolv_ctx->correction_solver, "Correction solver (with Schur factorization)"));
1280         PetscCall(PCShellSetContext(msolv_ctx->correction_solver, msolv_ctx));
1281         PetscCall(PCShellSetView(msolv_ctx->interior_solver, PCBDDCReuseSolvers_View));
1282         PetscCall(PCShellSetApply(msolv_ctx->correction_solver, PCBDDCReuseSolvers_Correction));
1283         PetscCall(PCShellSetApplyTranspose(msolv_ctx->correction_solver, PCBDDCReuseSolvers_CorrectionTranspose));
1284 
1285         /* scatter and vecs for Schur complement solver */
1286         PetscCall(MatCreateVecs(S_all, &msolv_ctx->sol_B, &msolv_ctx->rhs_B));
1287         PetscCall(MatCreateVecs(sub_schurs->S, &vec1_B, NULL));
1288         if (!schur_has_vertices) {
1289           PetscCall(ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap, IS_GTOLM_DROP, is_A_all, &msolv_ctx->is_B));
1290           PetscCall(VecScatterCreate(vec1_B, msolv_ctx->is_B, msolv_ctx->sol_B, NULL, &msolv_ctx->correction_scatter_B));
1291           PetscCall(PetscObjectReference((PetscObject)is_A_all));
1292           msolv_ctx->is_R = is_A_all;
1293         } else {
1294           IS              is_B_all;
1295           const PetscInt *idxs;
1296           PetscInt        dual, n_v, n;
1297 
1298           PetscCall(ISGetLocalSize(sub_schurs->is_vertices, &n_v));
1299           dual = size_schur - n_v;
1300           PetscCall(ISGetLocalSize(is_A_all, &n));
1301           PetscCall(ISGetIndices(is_A_all, &idxs));
1302           PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)is_A_all), dual, idxs + n_I, PETSC_COPY_VALUES, &is_B_all));
1303           PetscCall(ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap, IS_GTOLM_DROP, is_B_all, &msolv_ctx->is_B));
1304           PetscCall(ISDestroy(&is_B_all));
1305           PetscCall(ISCreateStride(PetscObjectComm((PetscObject)is_A_all), dual, 0, 1, &is_B_all));
1306           PetscCall(VecScatterCreate(vec1_B, msolv_ctx->is_B, msolv_ctx->sol_B, is_B_all, &msolv_ctx->correction_scatter_B));
1307           PetscCall(ISDestroy(&is_B_all));
1308           PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)is_A_all), n - n_v, idxs, PETSC_COPY_VALUES, &msolv_ctx->is_R));
1309           PetscCall(ISRestoreIndices(is_A_all, &idxs));
1310         }
1311         PetscCall(ISGetLocalSize(msolv_ctx->is_R, &n_R));
1312         PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF, n_R, n_R, 0, NULL, &Afake));
1313         PetscCall(MatAssemblyBegin(Afake, MAT_FINAL_ASSEMBLY));
1314         PetscCall(MatAssemblyEnd(Afake, MAT_FINAL_ASSEMBLY));
1315         PetscCall(PCSetOperators(msolv_ctx->correction_solver, Afake, Afake));
1316         PetscCall(MatDestroy(&Afake));
1317         PetscCall(VecDestroy(&vec1_B));
1318       }
1319       /* communicate benign info to solver context */
1320       if (benign_n) {
1321         PetscScalar *array;
1322 
1323         msolv_ctx->benign_n             = benign_n;
1324         msolv_ctx->benign_zerodiag_subs = is_p_r;
1325         PetscCall(PetscMalloc1(benign_n, &msolv_ctx->benign_save_vals));
1326         msolv_ctx->benign_csAIB = cs_AIB_mat;
1327         PetscCall(MatCreateVecs(cs_AIB_mat, &msolv_ctx->benign_corr_work, NULL));
1328         PetscCall(VecGetArray(msolv_ctx->benign_corr_work, &array));
1329         PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, size_schur, array, &msolv_ctx->benign_dummy_schur_vec));
1330         PetscCall(VecRestoreArray(msolv_ctx->benign_corr_work, &array));
1331         msolv_ctx->benign_AIIm1ones = benign_AIIm1_ones_mat;
1332       }
1333     } else {
1334       if (sub_schurs->reuse_solver) PetscCall(PCBDDCReuseSolversReset(sub_schurs->reuse_solver));
1335       PetscCall(PetscFree(sub_schurs->reuse_solver));
1336     }
1337     PetscCall(MatDestroy(&A));
1338     PetscCall(ISDestroy(&is_A_all));
1339 
1340     /* Work arrays */
1341     PetscCall(PetscMalloc1(max_subset_size * max_subset_size, &work));
1342 
1343     /* S_Ej_all */
1344     cum = cum2 = 0;
1345     PetscCall(MatDenseGetArrayRead(S_all, &rS_data));
1346     PetscCall(MatSeqAIJGetArray(sub_schurs->S_Ej_all, &SEj_arr));
1347     if (sub_schurs->sum_S_Ej_inv_all) PetscCall(MatSeqAIJGetArray(sub_schurs->sum_S_Ej_inv_all, &SEjinv_arr));
1348     if (sub_schurs->gdsw) PetscCall(MatCreateSubMatrices(sub_schurs->A, sub_schurs->n_subs, sub_schurs->is_subs, sub_schurs->is_subs, MAT_INITIAL_MATRIX, &gdswA));
1349     for (i = 0; i < sub_schurs->n_subs; i++) {
1350       PetscInt j;
1351 
1352       /* get S_E (or K^i_EE for GDSW) */
1353       PetscCall(ISGetLocalSize(sub_schurs->is_subs[i], &subset_size));
1354       if (sub_schurs->gdsw) {
1355         Mat T;
1356 
1357         PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, subset_size, subset_size, work, &T));
1358         PetscCall(MatConvert(gdswA[i], MATDENSE, MAT_REUSE_MATRIX, &T));
1359         PetscCall(MatDestroy(&T));
1360       } else {
1361         if (S_lower_triangular) { /* I need to expand the upper triangular data (column oriented) */
1362           PetscInt k;
1363           for (k = 0; k < subset_size; k++) {
1364             for (j = k; j < subset_size; j++) {
1365               work[k * subset_size + j] = rS_data[cum2 + k * size_schur + j];
1366               work[j * subset_size + k] = PetscConj(rS_data[cum2 + k * size_schur + j]);
1367             }
1368           }
1369         } else { /* just copy to workspace */
1370           PetscInt k;
1371           for (k = 0; k < subset_size; k++) {
1372             for (j = 0; j < subset_size; j++) work[k * subset_size + j] = rS_data[cum2 + k * size_schur + j];
1373           }
1374         }
1375       }
1376       /* insert S_E values */
1377       if (sub_schurs->change) {
1378         Mat change_sub, SEj, T;
1379 
1380         /* change basis */
1381         PetscCall(KSPGetOperators(sub_schurs->change[i], &change_sub, NULL));
1382         PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, subset_size, subset_size, work, &SEj));
1383         if (!sub_schurs->change_with_qr) { /* currently there's no support for PtAP with P SeqAIJ */
1384           Mat T2;
1385           PetscCall(MatTransposeMatMult(change_sub, SEj, MAT_INITIAL_MATRIX, 1.0, &T2));
1386           PetscCall(MatMatMult(T2, change_sub, MAT_INITIAL_MATRIX, 1.0, &T));
1387           PetscCall(MatConvert(T, MATSEQDENSE, MAT_INPLACE_MATRIX, &T));
1388           PetscCall(MatDestroy(&T2));
1389         } else {
1390           PetscCall(MatPtAP(SEj, change_sub, MAT_INITIAL_MATRIX, 1.0, &T));
1391         }
1392         PetscCall(MatCopy(T, SEj, SAME_NONZERO_PATTERN));
1393         PetscCall(MatDestroy(&T));
1394         PetscCall(MatZeroRowsColumnsIS(SEj, sub_schurs->change_primal_sub[i], 1.0, NULL, NULL));
1395         PetscCall(MatDestroy(&SEj));
1396       }
1397       PetscCall(PetscArraycpy(SEj_arr, work, subset_size * subset_size));
1398       if (compute_Stilda) {
1399         if (deluxe) { /* if adaptivity is requested, invert S_E blocks */
1400           Mat                M;
1401           const PetscScalar *vals;
1402           PetscBool          isdense, isdensecuda;
1403 
1404           PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, subset_size, subset_size, work, &M));
1405           PetscCall(MatSetOption(M, MAT_SPD, sub_schurs->is_posdef));
1406           PetscCall(MatSetOption(M, MAT_HERMITIAN, sub_schurs->is_hermitian));
1407           if (!PetscBTLookup(sub_schurs->is_edge, i)) PetscCall(MatSetType(M, Stype));
1408           PetscCall(PetscObjectTypeCompare((PetscObject)M, MATSEQDENSE, &isdense));
1409           PetscCall(PetscObjectTypeCompare((PetscObject)M, MATSEQDENSECUDA, &isdensecuda));
1410           if (use_cholesky) {
1411             PetscCall(MatCholeskyFactor(M, NULL, NULL));
1412           } else {
1413             PetscCall(MatLUFactor(M, NULL, NULL, NULL));
1414           }
1415           if (isdense) {
1416             PetscCall(MatSeqDenseInvertFactors_Private(M));
1417 #if defined(PETSC_HAVE_CUDA)
1418           } else if (isdensecuda) {
1419             PetscCall(MatSeqDenseCUDAInvertFactors_Internal(M));
1420 #endif
1421           } else SETERRQ(PetscObjectComm((PetscObject)M), PETSC_ERR_SUP, "Not implemented for type %s", Stype);
1422           PetscCall(MatDenseGetArrayRead(M, &vals));
1423           PetscCall(PetscArraycpy(SEjinv_arr, vals, subset_size * subset_size));
1424           PetscCall(MatDenseRestoreArrayRead(M, &vals));
1425           PetscCall(MatDestroy(&M));
1426         } else if (scaling) { /* not using deluxe */
1427           Mat          SEj;
1428           Vec          D;
1429           PetscScalar *array;
1430 
1431           PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, subset_size, subset_size, work, &SEj));
1432           PetscCall(VecGetArray(Dall, &array));
1433           PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, subset_size, array + cum, &D));
1434           PetscCall(VecRestoreArray(Dall, &array));
1435           PetscCall(VecShift(D, -1.));
1436           PetscCall(MatDiagonalScale(SEj, D, D));
1437           PetscCall(MatDestroy(&SEj));
1438           PetscCall(VecDestroy(&D));
1439           PetscCall(PetscArraycpy(SEj_arr, work, subset_size * subset_size));
1440         }
1441       }
1442       cum += subset_size;
1443       cum2 += subset_size * (size_schur + 1);
1444       SEj_arr += subset_size * subset_size;
1445       if (SEjinv_arr) SEjinv_arr += subset_size * subset_size;
1446     }
1447     if (sub_schurs->gdsw) PetscCall(MatDestroySubMatrices(sub_schurs->n_subs, &gdswA));
1448     PetscCall(MatDenseRestoreArrayRead(S_all, &rS_data));
1449     PetscCall(MatSeqAIJRestoreArray(sub_schurs->S_Ej_all, &SEj_arr));
1450     if (sub_schurs->sum_S_Ej_inv_all) PetscCall(MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_inv_all, &SEjinv_arr));
1451     if (solver_S) PetscCall(MatFactorRestoreSchurComplement(F, &S_all, MAT_FACTOR_SCHUR_UNFACTORED));
1452 
1453       /* may prevent from unneeded copies, since MUMPS or MKL_Pardiso always use CPU memory
1454        however, preliminary tests indicate using GPUs is still faster in the solve phase */
1455 #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
1456     if (reuse_solvers) {
1457       Mat                  St;
1458       MatFactorSchurStatus st;
1459 
1460       flg = PETSC_FALSE;
1461       PetscCall(PetscOptionsGetBool(NULL, sub_schurs->prefix, "-sub_schurs_schur_pin_to_cpu", &flg, NULL));
1462       PetscCall(MatFactorGetSchurComplement(F, &St, &st));
1463       PetscCall(MatBindToCPU(St, flg));
1464       PetscCall(MatFactorRestoreSchurComplement(F, &St, st));
1465     }
1466 #endif
1467 
1468     schur_factor = NULL;
1469     if (compute_Stilda && size_active_schur) {
1470       if (sub_schurs->n_subs == 1 && size_schur == size_active_schur && deluxe) { /* we already computed the inverse */
1471         PetscCall(MatSeqAIJGetArrayWrite(sub_schurs->sum_S_Ej_tilda_all, &SEjinv_arr));
1472         PetscCall(PetscArraycpy(SEjinv_arr, work, size_schur * size_schur));
1473         PetscCall(MatSeqAIJRestoreArrayWrite(sub_schurs->sum_S_Ej_tilda_all, &SEjinv_arr));
1474       } else {
1475         Mat S_all_inv = NULL;
1476 
1477         if (solver_S && !sub_schurs->gdsw) {
1478           /* for adaptive selection we need S^-1; for solver reusage we need S_\Delta\Delta^-1.
1479              The latter is not the principal subminor for S^-1. However, the factors can be reused since S_\Delta\Delta is the leading principal submatrix of S */
1480           if (factor_workaround) { /* invert without calling MatFactorInvertSchurComplement, since we are hacking */
1481             PetscScalar *data;
1482             PetscInt     nd = 0;
1483 
1484             PetscCheck(use_potr, PETSC_COMM_SELF, PETSC_ERR_SUP, "Factor update not yet implemented for non SPD matrices");
1485             PetscCall(MatFactorGetSchurComplement(F, &S_all_inv, NULL));
1486             PetscCall(MatDenseGetArray(S_all_inv, &data));
1487             if (sub_schurs->is_dir) { /* dirichlet dofs could have different scalings */
1488               PetscCall(ISGetLocalSize(sub_schurs->is_dir, &nd));
1489             }
1490 
1491             /* factor and invert activedofs and vertices (dirichlet dofs does not contribute) */
1492             if (schur_has_vertices) {
1493               Mat          M;
1494               PetscScalar *tdata;
1495               PetscInt     nv = 0, news;
1496 
1497               PetscCall(ISGetLocalSize(sub_schurs->is_vertices, &nv));
1498               news = size_active_schur + nv;
1499               PetscCall(PetscCalloc1(news * news, &tdata));
1500               for (i = 0; i < size_active_schur; i++) {
1501                 PetscCall(PetscArraycpy(tdata + i * (news + 1), data + i * (size_schur + 1), size_active_schur - i));
1502                 PetscCall(PetscArraycpy(tdata + i * (news + 1) + size_active_schur - i, data + i * size_schur + size_active_schur + nd, nv));
1503               }
1504               for (i = 0; i < nv; i++) {
1505                 PetscInt k = i + size_active_schur;
1506                 PetscCall(PetscArraycpy(tdata + k * (news + 1), data + (k + nd) * (size_schur + 1), nv - i));
1507               }
1508 
1509               PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, news, news, tdata, &M));
1510               PetscCall(MatSetOption(M, MAT_SPD, PETSC_TRUE));
1511               PetscCall(MatCholeskyFactor(M, NULL, NULL));
1512               /* save the factors */
1513               cum = 0;
1514               PetscCall(PetscMalloc1((size_active_schur * (size_active_schur + 1)) / 2 + nd, &schur_factor));
1515               for (i = 0; i < size_active_schur; i++) {
1516                 PetscCall(PetscArraycpy(schur_factor + cum, tdata + i * (news + 1), size_active_schur - i));
1517                 cum += size_active_schur - i;
1518               }
1519               for (i = 0; i < nd; i++) schur_factor[cum + i] = PetscSqrtReal(PetscRealPart(data[(i + size_active_schur) * (size_schur + 1)]));
1520               PetscCall(MatSeqDenseInvertFactors_Private(M));
1521               /* move back just the active dofs to the Schur complement */
1522               for (i = 0; i < size_active_schur; i++) PetscCall(PetscArraycpy(data + i * size_schur, tdata + i * news, size_active_schur));
1523               PetscCall(PetscFree(tdata));
1524               PetscCall(MatDestroy(&M));
1525             } else { /* we can factorize and invert just the activedofs */
1526               Mat          M;
1527               PetscScalar *aux;
1528 
1529               PetscCall(PetscMalloc1(nd, &aux));
1530               for (i = 0; i < nd; i++) aux[i] = 1.0 / data[(i + size_active_schur) * (size_schur + 1)];
1531               PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, size_active_schur, size_active_schur, data, &M));
1532               PetscCall(MatDenseSetLDA(M, size_schur));
1533               PetscCall(MatSetOption(M, MAT_SPD, PETSC_TRUE));
1534               PetscCall(MatCholeskyFactor(M, NULL, NULL));
1535               PetscCall(MatSeqDenseInvertFactors_Private(M));
1536               PetscCall(MatDestroy(&M));
1537               PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, size_schur, nd, data + size_active_schur * size_schur, &M));
1538               PetscCall(MatZeroEntries(M));
1539               PetscCall(MatDestroy(&M));
1540               PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, nd, size_schur, data + size_active_schur, &M));
1541               PetscCall(MatDenseSetLDA(M, size_schur));
1542               PetscCall(MatZeroEntries(M));
1543               PetscCall(MatDestroy(&M));
1544               for (i = 0; i < nd; i++) data[(i + size_active_schur) * (size_schur + 1)] = aux[i];
1545               PetscCall(PetscFree(aux));
1546             }
1547             PetscCall(MatDenseRestoreArray(S_all_inv, &data));
1548           } else { /* use MatFactor calls to invert S */
1549             PetscCall(MatFactorInvertSchurComplement(F));
1550             PetscCall(MatFactorGetSchurComplement(F, &S_all_inv, NULL));
1551           }
1552         } else if (!sub_schurs->gdsw) { /* we need to invert explicitly since we are not using MatFactor for S */
1553           PetscCall(PetscObjectReference((PetscObject)S_all));
1554           S_all_inv = S_all;
1555           PetscCall(MatDenseGetArray(S_all_inv, &S_data));
1556           PetscCall(PetscBLASIntCast(size_schur, &B_N));
1557           PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
1558           if (use_potr) {
1559             PetscCallBLAS("LAPACKpotrf", LAPACKpotrf_("L", &B_N, S_data, &B_N, &B_ierr));
1560             PetscCheck(!B_ierr, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in POTRF Lapack routine %d", (int)B_ierr);
1561             PetscCallBLAS("LAPACKpotri", LAPACKpotri_("L", &B_N, S_data, &B_N, &B_ierr));
1562             PetscCheck(!B_ierr, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in POTRI Lapack routine %d", (int)B_ierr);
1563           } else if (use_sytr) {
1564             PetscCallBLAS("LAPACKsytrf", LAPACKsytrf_("L", &B_N, S_data, &B_N, pivots, Bwork, &B_lwork, &B_ierr));
1565             PetscCheck(!B_ierr, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in SYTRF Lapack routine %d", (int)B_ierr);
1566             PetscCallBLAS("LAPACKsytri", LAPACKsytri_("L", &B_N, S_data, &B_N, pivots, Bwork, &B_ierr));
1567             PetscCheck(!B_ierr, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in SYTRI Lapack routine %d", (int)B_ierr);
1568           } else {
1569             PetscCallBLAS("LAPACKgetrf", LAPACKgetrf_(&B_N, &B_N, S_data, &B_N, pivots, &B_ierr));
1570             PetscCheck(!B_ierr, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in GETRF Lapack routine %d", (int)B_ierr);
1571             PetscCallBLAS("LAPACKgetri", LAPACKgetri_(&B_N, S_data, &B_N, pivots, Bwork, &B_lwork, &B_ierr));
1572             PetscCheck(!B_ierr, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in GETRI Lapack routine %d", (int)B_ierr);
1573           }
1574           PetscCall(PetscLogFlops(1.0 * size_schur * size_schur * size_schur));
1575           PetscCall(PetscFPTrapPop());
1576           PetscCall(MatDenseRestoreArray(S_all_inv, &S_data));
1577         } else if (sub_schurs->gdsw) {
1578           Mat      tS, tX, SEj, S_II, S_IE, S_EE;
1579           KSP      pS_II;
1580           PC       pS_II_pc;
1581           IS       EE, II;
1582           PetscInt nS;
1583 
1584           PetscCall(MatFactorCreateSchurComplement(F, &tS, NULL));
1585           PetscCall(MatGetSize(tS, &nS, NULL));
1586           PetscCall(MatSeqAIJGetArray(sub_schurs->sum_S_Ej_tilda_all, &SEjinv_arr));
1587           for (i = 0, cum = 0; i < sub_schurs->n_subs; i++) { /* naive implementation */
1588             PetscCall(ISGetLocalSize(sub_schurs->is_subs[i], &subset_size));
1589             PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, subset_size, subset_size, SEjinv_arr, &SEj));
1590 
1591             PetscCall(ISCreateStride(PETSC_COMM_SELF, subset_size, cum, 1, &EE));
1592             PetscCall(ISComplement(EE, 0, nS, &II));
1593             PetscCall(MatCreateSubMatrix(tS, II, II, MAT_INITIAL_MATRIX, &S_II));
1594             PetscCall(MatCreateSubMatrix(tS, II, EE, MAT_INITIAL_MATRIX, &S_IE));
1595             PetscCall(MatCreateSubMatrix(tS, EE, EE, MAT_INITIAL_MATRIX, &S_EE));
1596             PetscCall(ISDestroy(&II));
1597             PetscCall(ISDestroy(&EE));
1598 
1599             PetscCall(KSPCreate(PETSC_COMM_SELF, &pS_II));
1600             PetscCall(KSPSetType(pS_II, KSPPREONLY));
1601             PetscCall(KSPGetPC(pS_II, &pS_II_pc));
1602             PetscCall(PCSetType(pS_II_pc, PCSVD));
1603             PetscCall(KSPSetOptionsPrefix(pS_II, sub_schurs->prefix));
1604             PetscCall(KSPAppendOptionsPrefix(pS_II, "pseudo_"));
1605             PetscCall(KSPSetOperators(pS_II, S_II, S_II));
1606             PetscCall(MatDestroy(&S_II));
1607             PetscCall(KSPSetFromOptions(pS_II));
1608             PetscCall(KSPSetUp(pS_II));
1609             PetscCall(MatDuplicate(S_IE, MAT_DO_NOT_COPY_VALUES, &tX));
1610             PetscCall(KSPMatSolve(pS_II, S_IE, tX));
1611             PetscCall(KSPDestroy(&pS_II));
1612 
1613             PetscCall(MatTransposeMatMult(S_IE, tX, MAT_REUSE_MATRIX, PETSC_DEFAULT, &SEj));
1614             PetscCall(MatDestroy(&S_IE));
1615             PetscCall(MatDestroy(&tX));
1616             PetscCall(MatAYPX(SEj, -1, S_EE, SAME_NONZERO_PATTERN));
1617             PetscCall(MatDestroy(&S_EE));
1618 
1619             PetscCall(MatDestroy(&SEj));
1620             cum += subset_size;
1621             SEjinv_arr += subset_size * subset_size;
1622           }
1623           PetscCall(MatDestroy(&tS));
1624           PetscCall(MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_tilda_all, &SEjinv_arr));
1625         }
1626         /* S_Ej_tilda_all */
1627         cum = cum2 = 0;
1628         rS_data    = NULL;
1629         if (S_all_inv) PetscCall(MatDenseGetArrayRead(S_all_inv, &rS_data));
1630         PetscCall(MatSeqAIJGetArrayWrite(sub_schurs->sum_S_Ej_tilda_all, &SEjinv_arr));
1631         for (i = 0; i < sub_schurs->n_subs; i++) {
1632           PetscInt j;
1633 
1634           PetscCall(ISGetLocalSize(sub_schurs->is_subs[i], &subset_size));
1635           /* get (St^-1)_E */
1636           /* Unless we are changing the variables, I don't need to expand to upper triangular since St^-1
1637              will be properly accessed later during adaptive selection */
1638           if (rS_data) {
1639             if (S_lower_triangular) {
1640               PetscInt k;
1641               if (sub_schurs->change) {
1642                 for (k = 0; k < subset_size; k++) {
1643                   for (j = k; j < subset_size; j++) {
1644                     work[k * subset_size + j] = rS_data[cum2 + k * size_schur + j];
1645                     work[j * subset_size + k] = work[k * subset_size + j];
1646                   }
1647                 }
1648               } else {
1649                 for (k = 0; k < subset_size; k++) {
1650                   for (j = k; j < subset_size; j++) work[k * subset_size + j] = rS_data[cum2 + k * size_schur + j];
1651                 }
1652               }
1653             } else {
1654               PetscInt k;
1655               for (k = 0; k < subset_size; k++) {
1656                 for (j = 0; j < subset_size; j++) work[k * subset_size + j] = rS_data[cum2 + k * size_schur + j];
1657               }
1658             }
1659           }
1660           if (sub_schurs->change) {
1661             Mat         change_sub, SEj, T;
1662             PetscScalar val = sub_schurs->gdsw ? PETSC_SMALL : 1. / PETSC_SMALL;
1663 
1664             /* change basis */
1665             PetscCall(KSPGetOperators(sub_schurs->change[i], &change_sub, NULL));
1666             PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, subset_size, subset_size, rS_data ? work : SEjinv_arr, &SEj));
1667             if (!sub_schurs->change_with_qr) { /* currently there's no support for PtAP with P SeqAIJ */
1668               Mat T2;
1669               PetscCall(MatTransposeMatMult(change_sub, SEj, MAT_INITIAL_MATRIX, 1.0, &T2));
1670               PetscCall(MatMatMult(T2, change_sub, MAT_INITIAL_MATRIX, 1.0, &T));
1671               PetscCall(MatDestroy(&T2));
1672               PetscCall(MatConvert(T, MATSEQDENSE, MAT_INPLACE_MATRIX, &T));
1673             } else {
1674               PetscCall(MatPtAP(SEj, change_sub, MAT_INITIAL_MATRIX, 1.0, &T));
1675             }
1676             PetscCall(MatCopy(T, SEj, SAME_NONZERO_PATTERN));
1677             PetscCall(MatDestroy(&T));
1678             PetscCall(MatZeroRowsColumnsIS(SEj, sub_schurs->change_primal_sub[i], val, NULL, NULL));
1679             PetscCall(MatDestroy(&SEj));
1680           }
1681           if (rS_data) PetscCall(PetscArraycpy(SEjinv_arr, work, subset_size * subset_size));
1682           cum += subset_size;
1683           cum2 += subset_size * (size_schur + 1);
1684           SEjinv_arr += subset_size * subset_size;
1685         }
1686         PetscCall(MatSeqAIJRestoreArrayWrite(sub_schurs->sum_S_Ej_tilda_all, &SEjinv_arr));
1687         if (S_all_inv) {
1688           PetscCall(MatDenseRestoreArrayRead(S_all_inv, &rS_data));
1689           if (solver_S) {
1690             if (schur_has_vertices) {
1691               PetscCall(MatFactorRestoreSchurComplement(F, &S_all_inv, MAT_FACTOR_SCHUR_FACTORED));
1692             } else {
1693               PetscCall(MatFactorRestoreSchurComplement(F, &S_all_inv, MAT_FACTOR_SCHUR_INVERTED));
1694             }
1695           }
1696         }
1697         PetscCall(MatDestroy(&S_all_inv));
1698       }
1699 
1700       /* move back factors if needed */
1701       if (schur_has_vertices && factor_workaround && !sub_schurs->gdsw) {
1702         Mat          S_tmp;
1703         PetscInt     nd = 0;
1704         PetscScalar *data;
1705 
1706         PetscCheck(use_potr, PETSC_COMM_SELF, PETSC_ERR_SUP, "Factor update not yet implemented for non SPD matrices");
1707         PetscCheck(solver_S, PETSC_COMM_SELF, PETSC_ERR_PLIB, "This should not happen");
1708         PetscCall(MatFactorGetSchurComplement(F, &S_tmp, NULL));
1709         PetscCall(MatDenseGetArray(S_tmp, &data));
1710         PetscCall(PetscArrayzero(data, size_schur * size_schur));
1711 
1712         if (S_lower_triangular) {
1713           cum = 0;
1714           for (i = 0; i < size_active_schur; i++) {
1715             PetscCall(PetscArraycpy(data + i * (size_schur + 1), schur_factor + cum, size_active_schur - i));
1716             cum += size_active_schur - i;
1717           }
1718         } else {
1719           PetscCall(PetscArraycpy(data, schur_factor, size_schur * size_schur));
1720         }
1721         if (sub_schurs->is_dir) {
1722           PetscCall(ISGetLocalSize(sub_schurs->is_dir, &nd));
1723           for (i = 0; i < nd; i++) data[(i + size_active_schur) * (size_schur + 1)] = schur_factor[cum + i];
1724         }
1725         /* workaround: since I cannot modify the matrices used inside the solvers for the forward and backward substitutions,
1726              set the diagonal entry of the Schur factor to a very large value */
1727         for (i = size_active_schur + nd; i < size_schur; i++) data[i * (size_schur + 1)] = infty;
1728         PetscCall(MatDenseRestoreArray(S_tmp, &data));
1729         PetscCall(MatFactorRestoreSchurComplement(F, &S_tmp, MAT_FACTOR_SCHUR_FACTORED));
1730       }
1731     } else if (factor_workaround && !sub_schurs->gdsw) { /* we need to eliminate any unneeded coupling */
1732       PetscScalar *data;
1733       PetscInt     nd = 0;
1734 
1735       if (sub_schurs->is_dir) { /* dirichlet dofs could have different scalings */
1736         PetscCall(ISGetLocalSize(sub_schurs->is_dir, &nd));
1737       }
1738       PetscCall(MatFactorGetSchurComplement(F, &S_all, NULL));
1739       PetscCall(MatDenseGetArray(S_all, &data));
1740       for (i = 0; i < size_active_schur; i++) PetscCall(PetscArrayzero(data + i * size_schur + size_active_schur, size_schur - size_active_schur));
1741       for (i = size_active_schur + nd; i < size_schur; i++) {
1742         PetscCall(PetscArrayzero(data + i * size_schur + size_active_schur, size_schur - size_active_schur));
1743         data[i * (size_schur + 1)] = infty;
1744       }
1745       PetscCall(MatDenseRestoreArray(S_all, &data));
1746       PetscCall(MatFactorRestoreSchurComplement(F, &S_all, MAT_FACTOR_SCHUR_UNFACTORED));
1747     }
1748     PetscCall(PetscFree(work));
1749     PetscCall(PetscFree(schur_factor));
1750     PetscCall(VecDestroy(&Dall));
1751   }
1752   PetscCall(ISDestroy(&is_I_layer));
1753   PetscCall(MatDestroy(&S_all));
1754   PetscCall(MatDestroy(&A_BB));
1755   PetscCall(MatDestroy(&A_IB));
1756   PetscCall(MatDestroy(&A_BI));
1757   PetscCall(MatDestroy(&F));
1758 
1759   PetscCall(PetscMalloc1(sub_schurs->n_subs, &nnz));
1760   for (i = 0; i < sub_schurs->n_subs; i++) PetscCall(ISGetLocalSize(sub_schurs->is_subs[i], &nnz[i]));
1761   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, sub_schurs->n_subs, nnz, PETSC_OWN_POINTER, &is_I_layer));
1762   PetscCall(MatSetVariableBlockSizes(sub_schurs->S_Ej_all, sub_schurs->n_subs, nnz));
1763   PetscCall(MatAssemblyBegin(sub_schurs->S_Ej_all, MAT_FINAL_ASSEMBLY));
1764   PetscCall(MatAssemblyEnd(sub_schurs->S_Ej_all, MAT_FINAL_ASSEMBLY));
1765   if (compute_Stilda) {
1766     PetscCall(MatSetVariableBlockSizes(sub_schurs->sum_S_Ej_tilda_all, sub_schurs->n_subs, nnz));
1767     PetscCall(MatAssemblyBegin(sub_schurs->sum_S_Ej_tilda_all, MAT_FINAL_ASSEMBLY));
1768     PetscCall(MatAssemblyEnd(sub_schurs->sum_S_Ej_tilda_all, MAT_FINAL_ASSEMBLY));
1769     if (deluxe) {
1770       PetscCall(MatSetVariableBlockSizes(sub_schurs->sum_S_Ej_inv_all, sub_schurs->n_subs, nnz));
1771       PetscCall(MatAssemblyBegin(sub_schurs->sum_S_Ej_inv_all, MAT_FINAL_ASSEMBLY));
1772       PetscCall(MatAssemblyEnd(sub_schurs->sum_S_Ej_inv_all, MAT_FINAL_ASSEMBLY));
1773     }
1774   }
1775   PetscCall(ISDestroy(&is_I_layer));
1776 
1777   /* Get local part of (\sum_j S_Ej) */
1778   if (!sub_schurs->sum_S_Ej_all) PetscCall(MatDuplicate(sub_schurs->S_Ej_all, MAT_DO_NOT_COPY_VALUES, &sub_schurs->sum_S_Ej_all));
1779   PetscCall(VecSet(gstash, 0.0));
1780   PetscCall(MatSeqAIJGetArray(sub_schurs->S_Ej_all, &stasharray));
1781   PetscCall(VecPlaceArray(lstash, stasharray));
1782   PetscCall(VecScatterBegin(sstash, lstash, gstash, ADD_VALUES, SCATTER_FORWARD));
1783   PetscCall(VecScatterEnd(sstash, lstash, gstash, ADD_VALUES, SCATTER_FORWARD));
1784   PetscCall(MatSeqAIJRestoreArray(sub_schurs->S_Ej_all, &stasharray));
1785   PetscCall(VecResetArray(lstash));
1786   PetscCall(MatSeqAIJGetArray(sub_schurs->sum_S_Ej_all, &stasharray));
1787   PetscCall(VecPlaceArray(lstash, stasharray));
1788   PetscCall(VecScatterBegin(sstash, gstash, lstash, INSERT_VALUES, SCATTER_REVERSE));
1789   PetscCall(VecScatterEnd(sstash, gstash, lstash, INSERT_VALUES, SCATTER_REVERSE));
1790   PetscCall(MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_all, &stasharray));
1791   PetscCall(VecResetArray(lstash));
1792 
1793   /* Get local part of (\sum_j S^-1_Ej) (\sum_j St^-1_Ej) */
1794   if (compute_Stilda) {
1795     PetscCall(VecSet(gstash, 0.0));
1796     PetscCall(MatSeqAIJGetArray(sub_schurs->sum_S_Ej_tilda_all, &stasharray));
1797     PetscCall(VecPlaceArray(lstash, stasharray));
1798     PetscCall(VecScatterBegin(sstash, lstash, gstash, ADD_VALUES, SCATTER_FORWARD));
1799     PetscCall(VecScatterEnd(sstash, lstash, gstash, ADD_VALUES, SCATTER_FORWARD));
1800     PetscCall(VecScatterBegin(sstash, gstash, lstash, INSERT_VALUES, SCATTER_REVERSE));
1801     PetscCall(VecScatterEnd(sstash, gstash, lstash, INSERT_VALUES, SCATTER_REVERSE));
1802     PetscCall(MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_tilda_all, &stasharray));
1803     PetscCall(VecResetArray(lstash));
1804     if (deluxe) {
1805       PetscCall(VecSet(gstash, 0.0));
1806       PetscCall(MatSeqAIJGetArray(sub_schurs->sum_S_Ej_inv_all, &stasharray));
1807       PetscCall(VecPlaceArray(lstash, stasharray));
1808       PetscCall(VecScatterBegin(sstash, lstash, gstash, ADD_VALUES, SCATTER_FORWARD));
1809       PetscCall(VecScatterEnd(sstash, lstash, gstash, ADD_VALUES, SCATTER_FORWARD));
1810       PetscCall(VecScatterBegin(sstash, gstash, lstash, INSERT_VALUES, SCATTER_REVERSE));
1811       PetscCall(VecScatterEnd(sstash, gstash, lstash, INSERT_VALUES, SCATTER_REVERSE));
1812       PetscCall(MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_inv_all, &stasharray));
1813       PetscCall(VecResetArray(lstash));
1814     } else if (!sub_schurs->gdsw) {
1815       PetscScalar *array;
1816       PetscInt     cum;
1817 
1818       PetscCall(MatSeqAIJGetArray(sub_schurs->sum_S_Ej_tilda_all, &array));
1819       cum = 0;
1820       for (i = 0; i < sub_schurs->n_subs; i++) {
1821         PetscCall(ISGetLocalSize(sub_schurs->is_subs[i], &subset_size));
1822         PetscCall(PetscBLASIntCast(subset_size, &B_N));
1823         PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
1824         if (use_potr) {
1825           PetscCallBLAS("LAPACKpotrf", LAPACKpotrf_("L", &B_N, array + cum, &B_N, &B_ierr));
1826           PetscCheck(!B_ierr, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in POTRF Lapack routine %d", (int)B_ierr);
1827           PetscCallBLAS("LAPACKpotri", LAPACKpotri_("L", &B_N, array + cum, &B_N, &B_ierr));
1828           PetscCheck(!B_ierr, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in POTRI Lapack routine %d", (int)B_ierr);
1829         } else if (use_sytr) {
1830           PetscCallBLAS("LAPACKsytrf", LAPACKsytrf_("L", &B_N, array + cum, &B_N, pivots, Bwork, &B_lwork, &B_ierr));
1831           PetscCheck(!B_ierr, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in SYTRF Lapack routine %d", (int)B_ierr);
1832           PetscCallBLAS("LAPACKsytri", LAPACKsytri_("L", &B_N, array + cum, &B_N, pivots, Bwork, &B_ierr));
1833           PetscCheck(!B_ierr, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in SYTRI Lapack routine %d", (int)B_ierr);
1834         } else {
1835           PetscCallBLAS("LAPACKgetrf", LAPACKgetrf_(&B_N, &B_N, array + cum, &B_N, pivots, &B_ierr));
1836           PetscCheck(!B_ierr, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in GETRF Lapack routine %d", (int)B_ierr);
1837           PetscCallBLAS("LAPACKgetri", LAPACKgetri_(&B_N, array + cum, &B_N, pivots, Bwork, &B_lwork, &B_ierr));
1838           PetscCheck(!B_ierr, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in GETRI Lapack routine %d", (int)B_ierr);
1839         }
1840         PetscCall(PetscLogFlops(1.0 * subset_size * subset_size * subset_size));
1841         PetscCall(PetscFPTrapPop());
1842         cum += subset_size * subset_size;
1843       }
1844       PetscCall(MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_tilda_all, &array));
1845       PetscCall(PetscObjectReference((PetscObject)sub_schurs->sum_S_Ej_all));
1846       PetscCall(MatDestroy(&sub_schurs->sum_S_Ej_inv_all));
1847       sub_schurs->sum_S_Ej_inv_all = sub_schurs->sum_S_Ej_all;
1848     }
1849   }
1850   PetscCall(VecDestroy(&lstash));
1851   PetscCall(VecDestroy(&gstash));
1852   PetscCall(VecScatterDestroy(&sstash));
1853 
1854   if (matl_dbg_viewer) {
1855     if (sub_schurs->S_Ej_all) {
1856       PetscCall(PetscObjectSetName((PetscObject)sub_schurs->S_Ej_all, "SE"));
1857       PetscCall(MatView(sub_schurs->S_Ej_all, matl_dbg_viewer));
1858     }
1859     if (sub_schurs->sum_S_Ej_all) {
1860       PetscCall(PetscObjectSetName((PetscObject)sub_schurs->sum_S_Ej_all, "SSE"));
1861       PetscCall(MatView(sub_schurs->sum_S_Ej_all, matl_dbg_viewer));
1862     }
1863     if (sub_schurs->sum_S_Ej_inv_all) {
1864       PetscCall(PetscObjectSetName((PetscObject)sub_schurs->sum_S_Ej_inv_all, "SSEm"));
1865       PetscCall(MatView(sub_schurs->sum_S_Ej_inv_all, matl_dbg_viewer));
1866     }
1867     if (sub_schurs->sum_S_Ej_tilda_all) {
1868       PetscCall(PetscObjectSetName((PetscObject)sub_schurs->sum_S_Ej_tilda_all, "SSEt"));
1869       PetscCall(MatView(sub_schurs->sum_S_Ej_tilda_all, matl_dbg_viewer));
1870     }
1871   }
1872 
1873   /* free workspace */
1874   if (matl_dbg_viewer) PetscCall(PetscViewerFlush(matl_dbg_viewer));
1875   if (sub_schurs->debug) PetscCallMPI(MPI_Barrier(comm_n));
1876   PetscCall(PetscViewerDestroy(&matl_dbg_viewer));
1877   PetscCall(PetscFree2(Bwork, pivots));
1878   PetscCall(PetscCommDestroy(&comm_n));
1879   PetscFunctionReturn(PETSC_SUCCESS);
1880 }
1881 
1882 PetscErrorCode PCBDDCSubSchursInit(PCBDDCSubSchurs sub_schurs, const char *prefix, IS is_I, IS is_B, PCBDDCGraph graph, ISLocalToGlobalMapping BtoNmap, PetscBool copycc, PetscBool gdsw)
1883 {
1884   IS       *faces, *edges, *all_cc, vertices;
1885   PetscInt  s, i, n_faces, n_edges, n_all_cc;
1886   PetscBool is_sorted, ispardiso, ismumps;
1887 
1888   PetscFunctionBegin;
1889   PetscCall(ISSorted(is_I, &is_sorted));
1890   PetscCheck(is_sorted, PetscObjectComm((PetscObject)is_I), PETSC_ERR_PLIB, "IS for I dofs should be shorted");
1891   PetscCall(ISSorted(is_B, &is_sorted));
1892   PetscCheck(is_sorted, PetscObjectComm((PetscObject)is_B), PETSC_ERR_PLIB, "IS for B dofs should be shorted");
1893 
1894   /* reset any previous data */
1895   PetscCall(PCBDDCSubSchursReset(sub_schurs));
1896 
1897   sub_schurs->gdsw = gdsw;
1898 
1899   /* get index sets for faces and edges (already sorted by global ordering) */
1900   PetscCall(PCBDDCGraphGetCandidatesIS(graph, &n_faces, &faces, &n_edges, &edges, &vertices));
1901   n_all_cc = n_faces + n_edges;
1902   PetscCall(PetscBTCreate(n_all_cc, &sub_schurs->is_edge));
1903   PetscCall(PetscMalloc1(n_all_cc, &all_cc));
1904   n_all_cc = 0;
1905   for (i = 0; i < n_faces; i++) {
1906     PetscCall(ISGetSize(faces[i], &s));
1907     if (!s) continue;
1908     if (copycc) {
1909       PetscCall(ISDuplicate(faces[i], &all_cc[n_all_cc]));
1910     } else {
1911       PetscCall(PetscObjectReference((PetscObject)faces[i]));
1912       all_cc[n_all_cc] = faces[i];
1913     }
1914     n_all_cc++;
1915   }
1916   for (i = 0; i < n_edges; i++) {
1917     PetscCall(ISGetSize(edges[i], &s));
1918     if (!s) continue;
1919     if (copycc) {
1920       PetscCall(ISDuplicate(edges[i], &all_cc[n_all_cc]));
1921     } else {
1922       PetscCall(PetscObjectReference((PetscObject)edges[i]));
1923       all_cc[n_all_cc] = edges[i];
1924     }
1925     PetscCall(PetscBTSet(sub_schurs->is_edge, n_all_cc));
1926     n_all_cc++;
1927   }
1928   PetscCall(PetscObjectReference((PetscObject)vertices));
1929   sub_schurs->is_vertices = vertices;
1930   PetscCall(PCBDDCGraphRestoreCandidatesIS(graph, &n_faces, &faces, &n_edges, &edges, &vertices));
1931   sub_schurs->is_dir = NULL;
1932   PetscCall(PCBDDCGraphGetDirichletDofsB(graph, &sub_schurs->is_dir));
1933 
1934   /* Determine if MatFactor can be used */
1935   PetscCall(PetscStrallocpy(prefix, &sub_schurs->prefix));
1936 #if defined(PETSC_HAVE_MUMPS)
1937   PetscCall(PetscStrncpy(sub_schurs->mat_solver_type, MATSOLVERMUMPS, sizeof(sub_schurs->mat_solver_type)));
1938 #elif defined(PETSC_HAVE_MKL_PARDISO)
1939   PetscCall(PetscStrncpy(sub_schurs->mat_solver_type, MATSOLVERMKL_PARDISO, sizeof(sub_schurs->mat_solver_type)));
1940 #else
1941   PetscCall(PetscStrncpy(sub_schurs->mat_solver_type, MATSOLVERPETSC, sizeof(sub_schurs->mat_solver_type)));
1942 #endif
1943 #if defined(PETSC_USE_COMPLEX)
1944   sub_schurs->is_hermitian = PETSC_FALSE; /* Hermitian Cholesky is not supported by PETSc and external packages */
1945 #else
1946   sub_schurs->is_hermitian = PETSC_TRUE;
1947 #endif
1948   sub_schurs->is_posdef     = PETSC_TRUE;
1949   sub_schurs->is_symmetric  = PETSC_TRUE;
1950   sub_schurs->debug         = PETSC_FALSE;
1951   sub_schurs->restrict_comm = PETSC_FALSE;
1952   PetscOptionsBegin(PetscObjectComm((PetscObject)graph->l2gmap), sub_schurs->prefix, "BDDC sub_schurs options", "PC");
1953   PetscCall(PetscOptionsString("-sub_schurs_mat_solver_type", "Specific direct solver to use", NULL, sub_schurs->mat_solver_type, sub_schurs->mat_solver_type, sizeof(sub_schurs->mat_solver_type), NULL));
1954   PetscCall(PetscOptionsBool("-sub_schurs_symmetric", "Symmetric problem", NULL, sub_schurs->is_symmetric, &sub_schurs->is_symmetric, NULL));
1955   PetscCall(PetscOptionsBool("-sub_schurs_hermitian", "Hermitian problem", NULL, sub_schurs->is_hermitian, &sub_schurs->is_hermitian, NULL));
1956   PetscCall(PetscOptionsBool("-sub_schurs_posdef", "Positive definite problem", NULL, sub_schurs->is_posdef, &sub_schurs->is_posdef, NULL));
1957   PetscCall(PetscOptionsBool("-sub_schurs_restrictcomm", "Restrict communicator on active processes only", NULL, sub_schurs->restrict_comm, &sub_schurs->restrict_comm, NULL));
1958   PetscCall(PetscOptionsBool("-sub_schurs_debug", "Debug output", NULL, sub_schurs->debug, &sub_schurs->debug, NULL));
1959   PetscOptionsEnd();
1960   PetscCall(PetscStrcmp(sub_schurs->mat_solver_type, MATSOLVERMUMPS, &ismumps));
1961   PetscCall(PetscStrcmp(sub_schurs->mat_solver_type, MATSOLVERMKL_PARDISO, &ispardiso));
1962   sub_schurs->schur_explicit = (PetscBool)(ispardiso || ismumps);
1963 
1964   /* for reals, symmetric and hermitian are synonims */
1965 #if !defined(PETSC_USE_COMPLEX)
1966   sub_schurs->is_symmetric = (PetscBool)(sub_schurs->is_symmetric && sub_schurs->is_hermitian);
1967   sub_schurs->is_hermitian = sub_schurs->is_symmetric;
1968 #endif
1969 
1970   PetscCall(PetscObjectReference((PetscObject)is_I));
1971   sub_schurs->is_I = is_I;
1972   PetscCall(PetscObjectReference((PetscObject)is_B));
1973   sub_schurs->is_B = is_B;
1974   PetscCall(PetscObjectReference((PetscObject)graph->l2gmap));
1975   sub_schurs->l2gmap = graph->l2gmap;
1976   PetscCall(PetscObjectReference((PetscObject)BtoNmap));
1977   sub_schurs->BtoNmap            = BtoNmap;
1978   sub_schurs->n_subs             = n_all_cc;
1979   sub_schurs->is_subs            = all_cc;
1980   sub_schurs->S_Ej_all           = NULL;
1981   sub_schurs->sum_S_Ej_all       = NULL;
1982   sub_schurs->sum_S_Ej_inv_all   = NULL;
1983   sub_schurs->sum_S_Ej_tilda_all = NULL;
1984   sub_schurs->is_Ej_all          = NULL;
1985   PetscFunctionReturn(PETSC_SUCCESS);
1986 }
1987 
1988 PetscErrorCode PCBDDCSubSchursCreate(PCBDDCSubSchurs *sub_schurs)
1989 {
1990   PCBDDCSubSchurs schurs_ctx;
1991 
1992   PetscFunctionBegin;
1993   PetscCall(PetscNew(&schurs_ctx));
1994   schurs_ctx->n_subs = 0;
1995   *sub_schurs        = schurs_ctx;
1996   PetscFunctionReturn(PETSC_SUCCESS);
1997 }
1998 
1999 PetscErrorCode PCBDDCSubSchursReset(PCBDDCSubSchurs sub_schurs)
2000 {
2001   PetscInt i;
2002 
2003   PetscFunctionBegin;
2004   if (!sub_schurs) PetscFunctionReturn(PETSC_SUCCESS);
2005   PetscCall(PetscFree(sub_schurs->prefix));
2006   PetscCall(MatDestroy(&sub_schurs->A));
2007   PetscCall(MatDestroy(&sub_schurs->S));
2008   PetscCall(ISDestroy(&sub_schurs->is_I));
2009   PetscCall(ISDestroy(&sub_schurs->is_B));
2010   PetscCall(ISLocalToGlobalMappingDestroy(&sub_schurs->l2gmap));
2011   PetscCall(ISLocalToGlobalMappingDestroy(&sub_schurs->BtoNmap));
2012   PetscCall(MatDestroy(&sub_schurs->S_Ej_all));
2013   PetscCall(MatDestroy(&sub_schurs->sum_S_Ej_all));
2014   PetscCall(MatDestroy(&sub_schurs->sum_S_Ej_inv_all));
2015   PetscCall(MatDestroy(&sub_schurs->sum_S_Ej_tilda_all));
2016   PetscCall(ISDestroy(&sub_schurs->is_Ej_all));
2017   PetscCall(ISDestroy(&sub_schurs->is_vertices));
2018   PetscCall(ISDestroy(&sub_schurs->is_dir));
2019   PetscCall(PetscBTDestroy(&sub_schurs->is_edge));
2020   for (i = 0; i < sub_schurs->n_subs; i++) PetscCall(ISDestroy(&sub_schurs->is_subs[i]));
2021   if (sub_schurs->n_subs) PetscCall(PetscFree(sub_schurs->is_subs));
2022   if (sub_schurs->reuse_solver) PetscCall(PCBDDCReuseSolversReset(sub_schurs->reuse_solver));
2023   PetscCall(PetscFree(sub_schurs->reuse_solver));
2024   if (sub_schurs->change) {
2025     for (i = 0; i < sub_schurs->n_subs; i++) {
2026       PetscCall(KSPDestroy(&sub_schurs->change[i]));
2027       PetscCall(ISDestroy(&sub_schurs->change_primal_sub[i]));
2028     }
2029   }
2030   PetscCall(PetscFree(sub_schurs->change));
2031   PetscCall(PetscFree(sub_schurs->change_primal_sub));
2032   sub_schurs->n_subs = 0;
2033   PetscFunctionReturn(PETSC_SUCCESS);
2034 }
2035 
2036 PetscErrorCode PCBDDCSubSchursDestroy(PCBDDCSubSchurs *sub_schurs)
2037 {
2038   PetscFunctionBegin;
2039   PetscCall(PCBDDCSubSchursReset(*sub_schurs));
2040   PetscCall(PetscFree(*sub_schurs));
2041   PetscFunctionReturn(PETSC_SUCCESS);
2042 }
2043 
2044 static inline PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt *queue_tip, PetscInt n_prev, PetscBT touched, PetscInt *xadj, PetscInt *adjncy, PetscInt *n_added)
2045 {
2046   PetscInt i, j, n;
2047 
2048   PetscFunctionBegin;
2049   n = 0;
2050   for (i = -n_prev; i < 0; i++) {
2051     PetscInt start_dof = queue_tip[i];
2052     for (j = xadj[start_dof]; j < xadj[start_dof + 1]; j++) {
2053       PetscInt dof = adjncy[j];
2054       if (!PetscBTLookup(touched, dof)) {
2055         PetscCall(PetscBTSet(touched, dof));
2056         queue_tip[n] = dof;
2057         n++;
2058       }
2059     }
2060   }
2061   *n_added = n;
2062   PetscFunctionReturn(PETSC_SUCCESS);
2063 }
2064