xref: /petsc/src/ksp/pc/impls/bddc/bddcfetidp.c (revision e91c04dfc8a52dee1965211bb1cc8e5bf775178f)
1 #include <petsc/private/pcbddcimpl.h>
2 #include <petsc/private/pcbddcprivateimpl.h>
3 #include <petscblaslapack.h>
4 
5 static PetscErrorCode MatMult_BDdelta_deluxe_nonred(Mat A, Vec x, Vec y)
6 {
7   BDdelta_DN ctx;
8 
9   PetscFunctionBegin;
10   PetscCall(MatShellGetContext(A, &ctx));
11   PetscCall(MatMultTranspose(ctx->BD, x, ctx->work));
12   PetscCall(KSPSolveTranspose(ctx->kBD, ctx->work, y));
13   /* No PC so cannot propagate up failure in KSPSolveTranspose() */
14   PetscFunctionReturn(PETSC_SUCCESS);
15 }
16 
17 static PetscErrorCode MatMultTranspose_BDdelta_deluxe_nonred(Mat A, Vec x, Vec y)
18 {
19   BDdelta_DN ctx;
20 
21   PetscFunctionBegin;
22   PetscCall(MatShellGetContext(A, &ctx));
23   PetscCall(KSPSolve(ctx->kBD, x, ctx->work));
24   /* No PC so cannot propagate up failure in KSPSolve() */
25   PetscCall(MatMult(ctx->BD, ctx->work, y));
26   PetscFunctionReturn(PETSC_SUCCESS);
27 }
28 
29 static PetscErrorCode MatDestroy_BDdelta_deluxe_nonred(Mat A)
30 {
31   BDdelta_DN ctx;
32 
33   PetscFunctionBegin;
34   PetscCall(MatShellGetContext(A, &ctx));
35   PetscCall(MatDestroy(&ctx->BD));
36   PetscCall(KSPDestroy(&ctx->kBD));
37   PetscCall(VecDestroy(&ctx->work));
38   PetscCall(PetscFree(ctx));
39   PetscFunctionReturn(PETSC_SUCCESS);
40 }
41 
42 PetscErrorCode PCBDDCCreateFETIDPMatContext(PC pc, FETIDPMat_ctx *fetidpmat_ctx)
43 {
44   FETIDPMat_ctx newctx;
45 
46   PetscFunctionBegin;
47   PetscCall(PetscNew(&newctx));
48   /* increase the reference count for BDDC preconditioner */
49   PetscCall(PetscObjectReference((PetscObject)pc));
50   newctx->pc     = pc;
51   *fetidpmat_ctx = newctx;
52   PetscFunctionReturn(PETSC_SUCCESS);
53 }
54 
55 PetscErrorCode PCBDDCCreateFETIDPPCContext(PC pc, FETIDPPC_ctx *fetidppc_ctx)
56 {
57   FETIDPPC_ctx newctx;
58 
59   PetscFunctionBegin;
60   PetscCall(PetscNew(&newctx));
61   /* increase the reference count for BDDC preconditioner */
62   PetscCall(PetscObjectReference((PetscObject)pc));
63   newctx->pc    = pc;
64   *fetidppc_ctx = newctx;
65   PetscFunctionReturn(PETSC_SUCCESS);
66 }
67 
68 PetscErrorCode PCBDDCDestroyFETIDPMat(Mat A)
69 {
70   FETIDPMat_ctx mat_ctx;
71 
72   PetscFunctionBegin;
73   PetscCall(MatShellGetContext(A, &mat_ctx));
74   PetscCall(VecDestroy(&mat_ctx->lambda_local));
75   PetscCall(VecDestroy(&mat_ctx->temp_solution_D));
76   PetscCall(VecDestroy(&mat_ctx->temp_solution_B));
77   PetscCall(MatDestroy(&mat_ctx->B_delta));
78   PetscCall(MatDestroy(&mat_ctx->B_Ddelta));
79   PetscCall(MatDestroy(&mat_ctx->B_BB));
80   PetscCall(MatDestroy(&mat_ctx->B_BI));
81   PetscCall(MatDestroy(&mat_ctx->Bt_BB));
82   PetscCall(MatDestroy(&mat_ctx->Bt_BI));
83   PetscCall(MatDestroy(&mat_ctx->C));
84   PetscCall(VecDestroy(&mat_ctx->rhs_flip));
85   PetscCall(VecDestroy(&mat_ctx->vP));
86   PetscCall(VecDestroy(&mat_ctx->xPg));
87   PetscCall(VecDestroy(&mat_ctx->yPg));
88   PetscCall(VecScatterDestroy(&mat_ctx->l2g_lambda));
89   PetscCall(VecScatterDestroy(&mat_ctx->l2g_lambda_only));
90   PetscCall(VecScatterDestroy(&mat_ctx->l2g_p));
91   PetscCall(VecScatterDestroy(&mat_ctx->g2g_p));
92   PetscCall(PCDestroy(&mat_ctx->pc)); /* decrease PCBDDC reference count */
93   PetscCall(ISDestroy(&mat_ctx->pressure));
94   PetscCall(ISDestroy(&mat_ctx->lagrange));
95   PetscCall(PetscFree(mat_ctx));
96   PetscFunctionReturn(PETSC_SUCCESS);
97 }
98 
99 PetscErrorCode PCBDDCDestroyFETIDPPC(PC pc)
100 {
101   FETIDPPC_ctx pc_ctx;
102 
103   PetscFunctionBegin;
104   PetscCall(PCShellGetContext(pc, &pc_ctx));
105   PetscCall(VecDestroy(&pc_ctx->lambda_local));
106   PetscCall(MatDestroy(&pc_ctx->B_Ddelta));
107   PetscCall(VecScatterDestroy(&pc_ctx->l2g_lambda));
108   PetscCall(MatDestroy(&pc_ctx->S_j));
109   PetscCall(PCDestroy(&pc_ctx->pc)); /* decrease PCBDDC reference count */
110   PetscCall(VecDestroy(&pc_ctx->xPg));
111   PetscCall(VecDestroy(&pc_ctx->yPg));
112   PetscCall(PetscFree(pc_ctx));
113   PetscFunctionReturn(PETSC_SUCCESS);
114 }
115 
116 PetscErrorCode PCBDDCSetupFETIDPMatContext(FETIDPMat_ctx fetidpmat_ctx)
117 {
118   PC_IS          *pcis      = (PC_IS *)fetidpmat_ctx->pc->data;
119   PC_BDDC        *pcbddc    = (PC_BDDC *)fetidpmat_ctx->pc->data;
120   PCBDDCGraph     mat_graph = pcbddc->mat_graph;
121   Mat_IS         *matis     = (Mat_IS *)fetidpmat_ctx->pc->pmat->data;
122   MPI_Comm        comm;
123   Mat             ScalingMat, BD1, BD2;
124   Vec             fetidp_global;
125   IS              IS_l2g_lambda;
126   IS              subset, subset_mult, subset_n, isvert;
127   PetscBool       skip_node, fully_redundant;
128   PetscInt        i, j, k, s, n_boundary_dofs, n_global_lambda, n_vertices, partial_sum;
129   PetscInt        cum, n_local_lambda, n_lambda_for_dof, dual_size, n_neg_values, n_pos_values, buf_size;
130   PetscMPIInt     rank, size;
131   PetscScalar     scalar_value;
132   const PetscInt *vertex_indices;
133   PetscInt       *dual_dofs_boundary_indices, *aux_local_numbering_1;
134   const PetscInt *aux_global_numbering;
135   PetscInt       *aux_sums, *cols_B_delta, *l2g_indices;
136   PetscScalar    *array, *scaling_factors, *vals_B_delta;
137   PetscScalar   **all_factors;
138   PetscInt       *aux_local_numbering_2;
139   PetscInt       *count, **neighbours_set;
140   PetscLayout     llay;
141 
142   /* saddlepoint */
143   ISLocalToGlobalMapping l2gmap_p;
144   PetscLayout            play;
145   IS                     gP, pP;
146   PetscInt               nPl, nPg, nPgl;
147 
148   PetscFunctionBegin;
149   PetscCall(PetscObjectGetComm((PetscObject)fetidpmat_ctx->pc, &comm));
150   PetscCallMPI(MPI_Comm_rank(comm, &rank));
151   PetscCallMPI(MPI_Comm_size(comm, &size));
152 
153   /* saddlepoint */
154   nPl      = 0;
155   nPg      = 0;
156   nPgl     = 0;
157   gP       = NULL;
158   pP       = NULL;
159   l2gmap_p = NULL;
160   play     = NULL;
161   PetscCall(PetscObjectQuery((PetscObject)fetidpmat_ctx->pc, "__KSPFETIDP_pP", (PetscObject *)&pP));
162   if (pP) { /* saddle point */
163     /* subdomain pressures in global numbering */
164     PetscCall(PetscObjectQuery((PetscObject)fetidpmat_ctx->pc, "__KSPFETIDP_gP", (PetscObject *)&gP));
165     PetscCheck(gP, PETSC_COMM_SELF, PETSC_ERR_PLIB, "gP not present");
166     PetscCall(ISGetLocalSize(gP, &nPl));
167     PetscCall(VecCreate(PETSC_COMM_SELF, &fetidpmat_ctx->vP));
168     PetscCall(VecSetSizes(fetidpmat_ctx->vP, nPl, nPl));
169     PetscCall(VecSetType(fetidpmat_ctx->vP, VECSTANDARD));
170     PetscCall(VecSetUp(fetidpmat_ctx->vP));
171 
172     /* pressure matrix */
173     PetscCall(PetscObjectQuery((PetscObject)fetidpmat_ctx->pc, "__KSPFETIDP_C", (PetscObject *)&fetidpmat_ctx->C));
174     if (!fetidpmat_ctx->C) { /* null pressure block, compute layout and global numbering for pressures */
175       IS Pg;
176 
177       PetscCall(ISRenumber(gP, NULL, &nPg, &Pg));
178       PetscCall(ISLocalToGlobalMappingCreateIS(Pg, &l2gmap_p));
179       PetscCall(ISDestroy(&Pg));
180       PetscCall(PetscLayoutCreate(comm, &play));
181       PetscCall(PetscLayoutSetBlockSize(play, 1));
182       PetscCall(PetscLayoutSetSize(play, nPg));
183       PetscCall(ISGetLocalSize(pP, &nPgl));
184       PetscCall(PetscLayoutSetLocalSize(play, nPgl));
185       PetscCall(PetscLayoutSetUp(play));
186     } else {
187       PetscCall(PetscObjectReference((PetscObject)fetidpmat_ctx->C));
188       PetscCall(MatISGetLocalToGlobalMapping(fetidpmat_ctx->C, &l2gmap_p, NULL));
189       PetscCall(PetscObjectReference((PetscObject)l2gmap_p));
190       PetscCall(MatGetSize(fetidpmat_ctx->C, &nPg, NULL));
191       PetscCall(MatGetLocalSize(fetidpmat_ctx->C, NULL, &nPgl));
192       PetscCall(MatGetLayouts(fetidpmat_ctx->C, NULL, &llay));
193       PetscCall(PetscLayoutReference(llay, &play));
194     }
195     PetscCall(VecCreateMPIWithArray(comm, 1, nPgl, nPg, NULL, &fetidpmat_ctx->xPg));
196     PetscCall(VecCreateMPIWithArray(comm, 1, nPgl, nPg, NULL, &fetidpmat_ctx->yPg));
197 
198     /* import matrices for pressures coupling */
199     PetscCall(PetscObjectQuery((PetscObject)fetidpmat_ctx->pc, "__KSPFETIDP_B_BI", (PetscObject *)&fetidpmat_ctx->B_BI));
200     PetscCheck(fetidpmat_ctx->B_BI, PETSC_COMM_SELF, PETSC_ERR_PLIB, "B_BI not present");
201     PetscCall(PetscObjectReference((PetscObject)fetidpmat_ctx->B_BI));
202 
203     PetscCall(PetscObjectQuery((PetscObject)fetidpmat_ctx->pc, "__KSPFETIDP_B_BB", (PetscObject *)&fetidpmat_ctx->B_BB));
204     PetscCheck(fetidpmat_ctx->B_BB, PETSC_COMM_SELF, PETSC_ERR_PLIB, "B_BB not present");
205     PetscCall(PetscObjectReference((PetscObject)fetidpmat_ctx->B_BB));
206 
207     PetscCall(PetscObjectQuery((PetscObject)fetidpmat_ctx->pc, "__KSPFETIDP_Bt_BI", (PetscObject *)&fetidpmat_ctx->Bt_BI));
208     PetscCheck(fetidpmat_ctx->Bt_BI, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Bt_BI not present");
209     PetscCall(PetscObjectReference((PetscObject)fetidpmat_ctx->Bt_BI));
210 
211     PetscCall(PetscObjectQuery((PetscObject)fetidpmat_ctx->pc, "__KSPFETIDP_Bt_BB", (PetscObject *)&fetidpmat_ctx->Bt_BB));
212     PetscCheck(fetidpmat_ctx->Bt_BB, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Bt_BB not present");
213     PetscCall(PetscObjectReference((PetscObject)fetidpmat_ctx->Bt_BB));
214 
215     PetscCall(PetscObjectQuery((PetscObject)fetidpmat_ctx->pc, "__KSPFETIDP_flip", (PetscObject *)&fetidpmat_ctx->rhs_flip));
216     if (fetidpmat_ctx->rhs_flip) PetscCall(PetscObjectReference((PetscObject)fetidpmat_ctx->rhs_flip));
217   }
218 
219   /* Default type of lagrange multipliers is non-redundant */
220   fully_redundant = fetidpmat_ctx->fully_redundant;
221 
222   /* Evaluate local and global number of lagrange multipliers */
223   PetscCall(VecSet(pcis->vec1_N, 0.0));
224   n_local_lambda  = 0;
225   partial_sum     = 0;
226   n_boundary_dofs = 0;
227 
228   /* Get Vertices used to define the BDDC */
229   PetscCall(PCBDDCGraphGetCandidatesIS(mat_graph, NULL, NULL, NULL, NULL, &isvert));
230   PetscCall(ISGetLocalSize(isvert, &n_vertices));
231   PetscCall(ISGetIndices(isvert, &vertex_indices));
232 
233   dual_size = pcis->n_B - n_vertices;
234   PetscCall(PetscMalloc1(dual_size, &dual_dofs_boundary_indices));
235   PetscCall(PetscMalloc1(dual_size, &aux_local_numbering_1));
236   PetscCall(PetscMalloc1(dual_size, &aux_local_numbering_2));
237 
238   /* the code below does not support multiple subdomains per process
239      error out in this case
240      TODO: I guess I can use PetscSFGetMultiSF and the code will be easier and more general */
241   PetscCall(PetscMalloc2(pcis->n, &count, pcis->n, &neighbours_set));
242   for (i = 0, j = 0; i < pcis->n; i++) j += mat_graph->nodes[i].count;
243   if (pcis->n) PetscCall(PetscMalloc1(j, &neighbours_set[0]));
244   for (i = 0; i < pcis->n; i++) {
245     PCBDDCGraphNode *node = &mat_graph->nodes[i];
246 
247     count[i] = 0;
248     for (j = 0; j < node->count; j++) {
249       if (node->neighbours_set[j] == rank) continue;
250       neighbours_set[i][count[i]++] = node->neighbours_set[j];
251     }
252     PetscCheck(count[i] == node->count - 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Multiple subdomains per process not supported");
253     s = count[i];
254     PetscCall(PetscSortRemoveDupsInt(count + i, neighbours_set[i]));
255     PetscCheck(s == count[i], PETSC_COMM_SELF, PETSC_ERR_SUP, "Multiple subdomains per process not supported");
256     if (i != pcis->n - 1) neighbours_set[i + 1] = neighbours_set[i] + count[i];
257   }
258 
259   PetscCall(VecGetArray(pcis->vec1_N, &array));
260   for (i = 0, s = 0; i < pcis->n; i++) {
261     j = count[i]; /* RECALL: count[i] does not count myself */
262     if (j > 0) n_boundary_dofs++;
263     skip_node = PETSC_FALSE;
264     if (s < n_vertices && vertex_indices[s] == i) { /* it works for a sorted set of vertices */
265       skip_node = PETSC_TRUE;
266       s++;
267     }
268     if (j < 1) skip_node = PETSC_TRUE;
269     if (mat_graph->nodes[i].special_dof == PCBDDCGRAPH_DIRICHLET_MARK) skip_node = PETSC_TRUE;
270     if (!skip_node) {
271       if (fully_redundant) {
272         /* fully redundant set of lagrange multipliers */
273         n_lambda_for_dof = (j * (j + 1)) / 2;
274       } else {
275         n_lambda_for_dof = j;
276       }
277       n_local_lambda += j;
278       /* needed to evaluate global number of lagrange multipliers */
279       array[i] = (1.0 * n_lambda_for_dof) / (j + 1.0); /* already scaled for the next global sum */
280       /* store some data needed */
281       dual_dofs_boundary_indices[partial_sum] = n_boundary_dofs - 1;
282       aux_local_numbering_1[partial_sum]      = i;
283       aux_local_numbering_2[partial_sum]      = n_lambda_for_dof;
284       partial_sum++;
285     }
286   }
287   PetscCall(VecRestoreArray(pcis->vec1_N, &array));
288   PetscCall(ISRestoreIndices(isvert, &vertex_indices));
289   PetscCall(PCBDDCGraphRestoreCandidatesIS(mat_graph, NULL, NULL, NULL, NULL, &isvert));
290   dual_size = partial_sum;
291 
292   /* compute global ordering of lagrange multipliers and associate l2g map */
293   PetscCall(ISCreateGeneral(comm, partial_sum, aux_local_numbering_1, PETSC_COPY_VALUES, &subset_n));
294   PetscCall(ISLocalToGlobalMappingApplyIS(pcis->mapping, subset_n, &subset));
295   PetscCall(ISDestroy(&subset_n));
296   PetscCall(ISCreateGeneral(comm, partial_sum, aux_local_numbering_2, PETSC_OWN_POINTER, &subset_mult));
297   PetscCall(ISRenumber(subset, subset_mult, &fetidpmat_ctx->n_lambda, &subset_n));
298   PetscCall(ISDestroy(&subset));
299 
300   if (PetscDefined(USE_DEBUG)) {
301     PetscCall(VecSet(pcis->vec1_global, 0.0));
302     PetscCall(VecScatterBegin(matis->rctx, pcis->vec1_N, pcis->vec1_global, ADD_VALUES, SCATTER_REVERSE));
303     PetscCall(VecScatterEnd(matis->rctx, pcis->vec1_N, pcis->vec1_global, ADD_VALUES, SCATTER_REVERSE));
304     PetscCall(VecSum(pcis->vec1_global, &scalar_value));
305     i = (PetscInt)PetscRealPart(scalar_value);
306     PetscCheck(i == fetidpmat_ctx->n_lambda, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Global number of multipliers mismatch! (%" PetscInt_FMT " != %" PetscInt_FMT ")", fetidpmat_ctx->n_lambda, i);
307   }
308 
309   /* init data for scaling factors exchange */
310   if (!pcbddc->use_deluxe_scaling) {
311     PetscInt    *ptrs_buffer, neigh_position;
312     PetscScalar *send_buffer, *recv_buffer;
313     MPI_Request *send_reqs, *recv_reqs;
314     PetscMPIInt  nreqs;
315 
316     partial_sum = 0;
317     PetscCall(PetscMalloc1(pcis->n_neigh, &ptrs_buffer));
318     PetscCall(PetscMalloc1(PetscMax(pcis->n_neigh - 1, 0), &send_reqs));
319     PetscCall(PetscMalloc1(PetscMax(pcis->n_neigh - 1, 0), &recv_reqs));
320     PetscCall(PetscMalloc1(pcis->n + 1, &all_factors));
321     if (pcis->n_neigh > 0) ptrs_buffer[0] = 0;
322     for (i = 1; i < pcis->n_neigh; i++) {
323       partial_sum += pcis->n_shared[i];
324       ptrs_buffer[i] = ptrs_buffer[i - 1] + pcis->n_shared[i];
325     }
326     PetscCall(PetscMalloc1(partial_sum, &send_buffer));
327     PetscCall(PetscMalloc1(partial_sum, &recv_buffer));
328     PetscCall(PetscMalloc1(partial_sum, &all_factors[0]));
329     for (i = 0; i < pcis->n - 1; i++) {
330       j                  = count[i];
331       all_factors[i + 1] = all_factors[i] + j;
332     }
333 
334     /* scatter B scaling to N vec */
335     PetscCall(VecScatterBegin(pcis->N_to_B, pcis->D, pcis->vec1_N, INSERT_VALUES, SCATTER_REVERSE));
336     PetscCall(VecScatterEnd(pcis->N_to_B, pcis->D, pcis->vec1_N, INSERT_VALUES, SCATTER_REVERSE));
337     /* communications */
338     PetscCall(VecGetArrayRead(pcis->vec1_N, (const PetscScalar **)&array));
339     for (i = 1; i < pcis->n_neigh; i++) {
340       for (j = 0; j < pcis->n_shared[i]; j++) send_buffer[ptrs_buffer[i - 1] + j] = array[pcis->shared[i][j]];
341       buf_size = ptrs_buffer[i] - ptrs_buffer[i - 1];
342       PetscCallMPI(MPIU_Isend(&send_buffer[ptrs_buffer[i - 1]], buf_size, MPIU_SCALAR, pcis->neigh[i], 0, comm, &send_reqs[i - 1]));
343       PetscCallMPI(MPIU_Irecv(&recv_buffer[ptrs_buffer[i - 1]], buf_size, MPIU_SCALAR, pcis->neigh[i], 0, comm, &recv_reqs[i - 1]));
344     }
345     PetscCall(VecRestoreArrayRead(pcis->vec1_N, (const PetscScalar **)&array));
346     PetscCall(PetscMPIIntCast(pcis->n_neigh - 1, &nreqs));
347     if (pcis->n_neigh > 0) PetscCallMPI(MPI_Waitall(nreqs, recv_reqs, MPI_STATUSES_IGNORE));
348     /* put values in correct places */
349     for (i = 1; i < pcis->n_neigh; i++) {
350       for (j = 0; j < pcis->n_shared[i]; j++) {
351         k              = pcis->shared[i][j];
352         neigh_position = 0;
353         while (neighbours_set[k][neigh_position] != pcis->neigh[i]) { neigh_position++; }
354         all_factors[k][neigh_position] = recv_buffer[ptrs_buffer[i - 1] + j];
355       }
356     }
357     if (pcis->n_neigh > 0) PetscCallMPI(MPI_Waitall(nreqs, send_reqs, MPI_STATUSES_IGNORE));
358     PetscCall(PetscFree(send_reqs));
359     PetscCall(PetscFree(recv_reqs));
360     PetscCall(PetscFree(send_buffer));
361     PetscCall(PetscFree(recv_buffer));
362     PetscCall(PetscFree(ptrs_buffer));
363   }
364 
365   /* Compute B and B_delta (local actions) */
366   PetscCall(PetscMalloc1(pcis->n_neigh, &aux_sums));
367   PetscCall(PetscMalloc1(n_local_lambda, &l2g_indices));
368   PetscCall(PetscMalloc1(n_local_lambda, &vals_B_delta));
369   PetscCall(PetscMalloc1(n_local_lambda, &cols_B_delta));
370   if (!pcbddc->use_deluxe_scaling) {
371     PetscCall(PetscMalloc1(n_local_lambda, &scaling_factors));
372   } else {
373     scaling_factors = NULL;
374     all_factors     = NULL;
375   }
376   PetscCall(ISGetIndices(subset_n, &aux_global_numbering));
377   partial_sum = 0;
378   cum         = 0;
379   for (i = 0; i < dual_size; i++) {
380     n_global_lambda = aux_global_numbering[cum];
381     j               = count[aux_local_numbering_1[i]];
382     aux_sums[0]     = 0;
383     for (s = 1; s < j; s++) aux_sums[s] = aux_sums[s - 1] + j - s + 1;
384     if (all_factors) array = all_factors[aux_local_numbering_1[i]];
385     n_neg_values = 0;
386     while (n_neg_values < j && neighbours_set[aux_local_numbering_1[i]][n_neg_values] < rank) { n_neg_values++; }
387     n_pos_values = j - n_neg_values;
388     if (fully_redundant) {
389       for (s = 0; s < n_neg_values; s++) {
390         l2g_indices[partial_sum + s]  = aux_sums[s] + n_neg_values - s - 1 + n_global_lambda;
391         cols_B_delta[partial_sum + s] = dual_dofs_boundary_indices[i];
392         vals_B_delta[partial_sum + s] = -1.0;
393         if (!pcbddc->use_deluxe_scaling) scaling_factors[partial_sum + s] = array[s];
394       }
395       for (s = 0; s < n_pos_values; s++) {
396         l2g_indices[partial_sum + s + n_neg_values]  = aux_sums[n_neg_values] + s + n_global_lambda;
397         cols_B_delta[partial_sum + s + n_neg_values] = dual_dofs_boundary_indices[i];
398         vals_B_delta[partial_sum + s + n_neg_values] = 1.0;
399         if (!pcbddc->use_deluxe_scaling) scaling_factors[partial_sum + s + n_neg_values] = array[s + n_neg_values];
400       }
401       partial_sum += j;
402     } else {
403       /* l2g_indices and default cols and vals of B_delta */
404       for (s = 0; s < j; s++) {
405         l2g_indices[partial_sum + s]  = n_global_lambda + s;
406         cols_B_delta[partial_sum + s] = dual_dofs_boundary_indices[i];
407         vals_B_delta[partial_sum + s] = 0.0;
408       }
409       /* B_delta */
410       if (n_neg_values > 0) { /* there's a rank next to me to the left */
411         vals_B_delta[partial_sum + n_neg_values - 1] = -1.0;
412       }
413       if (n_neg_values < j) { /* there's a rank next to me to the right */
414         vals_B_delta[partial_sum + n_neg_values] = 1.0;
415       }
416       /* scaling as in Klawonn-Widlund 1999 */
417       if (!pcbddc->use_deluxe_scaling) {
418         for (s = 0; s < n_neg_values; s++) {
419           scalar_value = 0.0;
420           for (k = 0; k < s + 1; k++) scalar_value += array[k];
421           scaling_factors[partial_sum + s] = -scalar_value;
422         }
423         for (s = 0; s < n_pos_values; s++) {
424           scalar_value = 0.0;
425           for (k = s + n_neg_values; k < j; k++) scalar_value += array[k];
426           scaling_factors[partial_sum + s + n_neg_values] = scalar_value;
427         }
428       }
429       partial_sum += j;
430     }
431     cum += aux_local_numbering_2[i];
432   }
433   PetscCall(ISRestoreIndices(subset_n, &aux_global_numbering));
434   PetscCall(ISDestroy(&subset_mult));
435   PetscCall(ISDestroy(&subset_n));
436   PetscCall(PetscFree(aux_sums));
437   PetscCall(PetscFree(aux_local_numbering_1));
438   PetscCall(PetscFree(dual_dofs_boundary_indices));
439   if (all_factors) {
440     PetscCall(PetscFree(all_factors[0]));
441     PetscCall(PetscFree(all_factors));
442   }
443   if (pcis->n) PetscCall(PetscFree(neighbours_set[0]));
444   PetscCall(PetscFree2(count, neighbours_set));
445 
446   /* Create local part of B_delta */
447   PetscCall(MatCreate(PETSC_COMM_SELF, &fetidpmat_ctx->B_delta));
448   PetscCall(MatSetSizes(fetidpmat_ctx->B_delta, n_local_lambda, pcis->n_B, n_local_lambda, pcis->n_B));
449   PetscCall(MatSetType(fetidpmat_ctx->B_delta, MATSEQAIJ));
450   PetscCall(MatSeqAIJSetPreallocation(fetidpmat_ctx->B_delta, 1, NULL));
451   PetscCall(MatSetOption(fetidpmat_ctx->B_delta, MAT_IGNORE_ZERO_ENTRIES, PETSC_TRUE));
452   for (i = 0; i < n_local_lambda; i++) PetscCall(MatSetValue(fetidpmat_ctx->B_delta, i, cols_B_delta[i], vals_B_delta[i], INSERT_VALUES));
453   PetscCall(PetscFree(vals_B_delta));
454   PetscCall(MatAssemblyBegin(fetidpmat_ctx->B_delta, MAT_FINAL_ASSEMBLY));
455   PetscCall(MatAssemblyEnd(fetidpmat_ctx->B_delta, MAT_FINAL_ASSEMBLY));
456 
457   BD1 = NULL;
458   BD2 = NULL;
459   if (fully_redundant) {
460     PetscCheck(!pcbddc->use_deluxe_scaling, comm, PETSC_ERR_SUP, "Deluxe FETIDP with fully-redundant multipliers to be implemented");
461     PetscCall(MatCreate(PETSC_COMM_SELF, &ScalingMat));
462     PetscCall(MatSetSizes(ScalingMat, n_local_lambda, n_local_lambda, n_local_lambda, n_local_lambda));
463     PetscCall(MatSetType(ScalingMat, MATSEQAIJ));
464     PetscCall(MatSeqAIJSetPreallocation(ScalingMat, 1, NULL));
465     for (i = 0; i < n_local_lambda; i++) PetscCall(MatSetValue(ScalingMat, i, i, scaling_factors[i], INSERT_VALUES));
466     PetscCall(MatAssemblyBegin(ScalingMat, MAT_FINAL_ASSEMBLY));
467     PetscCall(MatAssemblyEnd(ScalingMat, MAT_FINAL_ASSEMBLY));
468     PetscCall(MatMatMult(ScalingMat, fetidpmat_ctx->B_delta, MAT_INITIAL_MATRIX, PETSC_DETERMINE, &fetidpmat_ctx->B_Ddelta));
469     PetscCall(MatDestroy(&ScalingMat));
470   } else {
471     PetscCall(MatCreate(PETSC_COMM_SELF, &fetidpmat_ctx->B_Ddelta));
472     PetscCall(MatSetSizes(fetidpmat_ctx->B_Ddelta, n_local_lambda, pcis->n_B, n_local_lambda, pcis->n_B));
473     if (!pcbddc->use_deluxe_scaling || !pcbddc->sub_schurs) {
474       PetscCall(MatSetType(fetidpmat_ctx->B_Ddelta, MATSEQAIJ));
475       PetscCall(MatSeqAIJSetPreallocation(fetidpmat_ctx->B_Ddelta, 1, NULL));
476       for (i = 0; i < n_local_lambda; i++) PetscCall(MatSetValue(fetidpmat_ctx->B_Ddelta, i, cols_B_delta[i], scaling_factors[i], INSERT_VALUES));
477       PetscCall(MatAssemblyBegin(fetidpmat_ctx->B_Ddelta, MAT_FINAL_ASSEMBLY));
478       PetscCall(MatAssemblyEnd(fetidpmat_ctx->B_Ddelta, MAT_FINAL_ASSEMBLY));
479     } else {
480       /* scaling as in Klawonn-Widlund 1999 */
481       PCBDDCDeluxeScaling deluxe_ctx = pcbddc->deluxe_ctx;
482       PCBDDCSubSchurs     sub_schurs = pcbddc->sub_schurs;
483       Mat                 T;
484       PetscScalar        *W, lwork, *Bwork;
485       const PetscInt     *idxs = NULL;
486       PetscInt            cum, mss, *nnz;
487       PetscBLASInt       *pivots, B_lwork, B_N, B_ierr;
488 
489       PetscCheck(pcbddc->deluxe_singlemat, comm, PETSC_ERR_USER, "Cannot compute B_Ddelta! rerun with -pc_bddc_deluxe_singlemat");
490       mss = 0;
491       PetscCall(PetscCalloc1(pcis->n_B, &nnz));
492       if (sub_schurs->is_Ej_all) {
493         PetscCall(ISGetIndices(sub_schurs->is_Ej_all, &idxs));
494         for (i = 0, cum = 0; i < sub_schurs->n_subs; i++) {
495           PetscInt subset_size;
496 
497           PetscCall(ISGetLocalSize(sub_schurs->is_subs[i], &subset_size));
498           for (j = 0; j < subset_size; j++) nnz[idxs[j + cum]] = subset_size;
499           mss = PetscMax(mss, subset_size);
500           cum += subset_size;
501         }
502       }
503       PetscCall(MatCreate(PETSC_COMM_SELF, &T));
504       PetscCall(MatSetSizes(T, pcis->n_B, pcis->n_B, pcis->n_B, pcis->n_B));
505       PetscCall(MatSetType(T, MATSEQAIJ));
506       PetscCall(MatSeqAIJSetPreallocation(T, 0, nnz));
507       PetscCall(PetscFree(nnz));
508 
509       /* workspace allocation */
510       B_lwork = 0;
511       if (mss) {
512         PetscScalar dummy = 1;
513 
514         B_lwork = -1;
515         PetscCall(PetscBLASIntCast(mss, &B_N));
516         PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
517         PetscCallBLAS("LAPACKgetri", LAPACKgetri_(&B_N, &dummy, &B_N, &B_N, &lwork, &B_lwork, &B_ierr));
518         PetscCall(PetscFPTrapPop());
519         PetscCheck(!B_ierr, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in query to GETRI Lapack routine %" PetscBLASInt_FMT, B_ierr);
520         PetscCall(PetscBLASIntCast((PetscInt)PetscRealPart(lwork), &B_lwork));
521       }
522       PetscCall(PetscMalloc3(mss * mss, &W, mss, &pivots, B_lwork, &Bwork));
523 
524       for (i = 0, cum = 0; i < sub_schurs->n_subs; i++) {
525         const PetscScalar *M;
526         PetscInt           subset_size;
527 
528         PetscCall(ISGetLocalSize(sub_schurs->is_subs[i], &subset_size));
529         PetscCall(PetscBLASIntCast(subset_size, &B_N));
530         PetscCall(MatDenseGetArrayRead(deluxe_ctx->seq_mat[i], &M));
531         PetscCall(PetscArraycpy(W, M, subset_size * subset_size));
532         PetscCall(MatDenseRestoreArrayRead(deluxe_ctx->seq_mat[i], &M));
533         PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
534         PetscCallBLAS("LAPACKgetrf", LAPACKgetrf_(&B_N, &B_N, W, &B_N, pivots, &B_ierr));
535         PetscCheck(!B_ierr, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in GETRF Lapack routine %" PetscBLASInt_FMT, B_ierr);
536         PetscCallBLAS("LAPACKgetri", LAPACKgetri_(&B_N, W, &B_N, pivots, Bwork, &B_lwork, &B_ierr));
537         PetscCheck(!B_ierr, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in GETRI Lapack routine %" PetscBLASInt_FMT, B_ierr);
538         PetscCall(PetscFPTrapPop());
539         /* silent static analyzer */
540         PetscCheck(idxs, PETSC_COMM_SELF, PETSC_ERR_PLIB, "IDXS not present");
541         PetscCall(MatSetValues(T, subset_size, idxs + cum, subset_size, idxs + cum, W, INSERT_VALUES));
542         cum += subset_size;
543       }
544       PetscCall(MatAssemblyBegin(T, MAT_FINAL_ASSEMBLY));
545       PetscCall(MatAssemblyEnd(T, MAT_FINAL_ASSEMBLY));
546       PetscCall(MatMatTransposeMult(T, fetidpmat_ctx->B_delta, MAT_INITIAL_MATRIX, PETSC_DETERMINE, &BD1));
547       PetscCall(MatMatMult(fetidpmat_ctx->B_delta, BD1, MAT_INITIAL_MATRIX, PETSC_DETERMINE, &BD2));
548       PetscCall(MatDestroy(&T));
549       PetscCall(PetscFree3(W, pivots, Bwork));
550       if (sub_schurs->is_Ej_all) PetscCall(ISRestoreIndices(sub_schurs->is_Ej_all, &idxs));
551     }
552   }
553   PetscCall(PetscFree(scaling_factors));
554   PetscCall(PetscFree(cols_B_delta));
555 
556   /* Layout of multipliers */
557   PetscCall(PetscLayoutCreate(comm, &llay));
558   PetscCall(PetscLayoutSetBlockSize(llay, 1));
559   PetscCall(PetscLayoutSetSize(llay, fetidpmat_ctx->n_lambda));
560   PetscCall(PetscLayoutSetUp(llay));
561   PetscCall(PetscLayoutGetLocalSize(llay, &fetidpmat_ctx->n));
562 
563   /* Local work vector of multipliers */
564   PetscCall(VecCreate(PETSC_COMM_SELF, &fetidpmat_ctx->lambda_local));
565   PetscCall(VecSetSizes(fetidpmat_ctx->lambda_local, n_local_lambda, n_local_lambda));
566   PetscCall(VecSetType(fetidpmat_ctx->lambda_local, VECSEQ));
567 
568   if (BD2) {
569     ISLocalToGlobalMapping l2g;
570     Mat                    T, TA, *pT;
571     IS                     is;
572     PetscInt               nl, N;
573     BDdelta_DN             ctx;
574 
575     PetscCall(PetscLayoutGetLocalSize(llay, &nl));
576     PetscCall(PetscLayoutGetSize(llay, &N));
577     PetscCall(MatCreate(comm, &T));
578     PetscCall(MatSetSizes(T, nl, nl, N, N));
579     PetscCall(MatSetType(T, MATIS));
580     PetscCall(ISLocalToGlobalMappingCreate(comm, 1, n_local_lambda, l2g_indices, PETSC_COPY_VALUES, &l2g));
581     PetscCall(MatSetLocalToGlobalMapping(T, l2g, l2g));
582     PetscCall(ISLocalToGlobalMappingDestroy(&l2g));
583     PetscCall(MatISSetLocalMat(T, BD2));
584     PetscCall(MatAssemblyBegin(T, MAT_FINAL_ASSEMBLY));
585     PetscCall(MatAssemblyEnd(T, MAT_FINAL_ASSEMBLY));
586     PetscCall(MatDestroy(&BD2));
587     PetscCall(MatConvert(T, MATAIJ, MAT_INITIAL_MATRIX, &TA));
588     PetscCall(MatDestroy(&T));
589     PetscCall(ISCreateGeneral(comm, n_local_lambda, l2g_indices, PETSC_USE_POINTER, &is));
590     PetscCall(MatCreateSubMatrices(TA, 1, &is, &is, MAT_INITIAL_MATRIX, &pT));
591     PetscCall(MatDestroy(&TA));
592     PetscCall(ISDestroy(&is));
593     BD2 = pT[0];
594     PetscCall(PetscFree(pT));
595 
596     /* B_Ddelta for non-redundant multipliers with deluxe scaling */
597     PetscCall(PetscNew(&ctx));
598     PetscCall(MatSetType(fetidpmat_ctx->B_Ddelta, MATSHELL));
599     PetscCall(MatShellSetContext(fetidpmat_ctx->B_Ddelta, ctx));
600     PetscCall(MatShellSetOperation(fetidpmat_ctx->B_Ddelta, MATOP_MULT, (void (*)(void))MatMult_BDdelta_deluxe_nonred));
601     PetscCall(MatShellSetOperation(fetidpmat_ctx->B_Ddelta, MATOP_MULT_TRANSPOSE, (void (*)(void))MatMultTranspose_BDdelta_deluxe_nonred));
602     PetscCall(MatShellSetOperation(fetidpmat_ctx->B_Ddelta, MATOP_DESTROY, (void (*)(void))MatDestroy_BDdelta_deluxe_nonred));
603     PetscCall(MatSetUp(fetidpmat_ctx->B_Ddelta));
604 
605     PetscCall(PetscObjectReference((PetscObject)BD1));
606     ctx->BD = BD1;
607     PetscCall(KSPCreate(PETSC_COMM_SELF, &ctx->kBD));
608     PetscCall(KSPSetNestLevel(ctx->kBD, fetidpmat_ctx->pc->kspnestlevel));
609     PetscCall(KSPSetOperators(ctx->kBD, BD2, BD2));
610     PetscCall(VecDuplicate(fetidpmat_ctx->lambda_local, &ctx->work));
611     fetidpmat_ctx->deluxe_nonred = PETSC_TRUE;
612   }
613   PetscCall(MatDestroy(&BD1));
614   PetscCall(MatDestroy(&BD2));
615 
616   /* fetidpmat sizes */
617   fetidpmat_ctx->n += nPgl;
618   fetidpmat_ctx->N = fetidpmat_ctx->n_lambda + nPg;
619 
620   /* Global vector for FETI-DP linear system */
621   PetscCall(VecCreate(comm, &fetidp_global));
622   PetscCall(VecSetSizes(fetidp_global, fetidpmat_ctx->n, fetidpmat_ctx->N));
623   PetscCall(VecSetType(fetidp_global, VECMPI));
624   PetscCall(VecSetUp(fetidp_global));
625 
626   /* Decide layout for fetidp dofs: if it is a saddle point problem
627      pressure is ordered first in the local part of the global vector
628      of the FETI-DP linear system */
629   if (nPg) {
630     Vec             v;
631     IS              IS_l2g_p, ais;
632     PetscLayout     alay;
633     const PetscInt *idxs, *pranges, *aranges, *lranges;
634     PetscInt       *l2g_indices_p, rst;
635     PetscMPIInt     rank;
636 
637     PetscCall(PetscMalloc1(nPl, &l2g_indices_p));
638     PetscCall(VecGetLayout(fetidp_global, &alay));
639     PetscCall(PetscLayoutGetRanges(alay, &aranges));
640     PetscCall(PetscLayoutGetRanges(play, &pranges));
641     PetscCall(PetscLayoutGetRanges(llay, &lranges));
642 
643     PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)fetidp_global), &rank));
644     PetscCall(ISCreateStride(PetscObjectComm((PetscObject)fetidp_global), pranges[rank + 1] - pranges[rank], aranges[rank], 1, &fetidpmat_ctx->pressure));
645     PetscCall(PetscObjectSetName((PetscObject)fetidpmat_ctx->pressure, "F_P"));
646     PetscCall(ISCreateStride(PetscObjectComm((PetscObject)fetidp_global), lranges[rank + 1] - lranges[rank], aranges[rank] + pranges[rank + 1] - pranges[rank], 1, &fetidpmat_ctx->lagrange));
647     PetscCall(PetscObjectSetName((PetscObject)fetidpmat_ctx->lagrange, "F_L"));
648     PetscCall(ISLocalToGlobalMappingGetIndices(l2gmap_p, &idxs));
649     /* shift local to global indices for pressure */
650     for (i = 0; i < nPl; i++) {
651       PetscMPIInt owner;
652 
653       PetscCall(PetscLayoutFindOwner(play, idxs[i], &owner));
654       l2g_indices_p[i] = idxs[i] - pranges[owner] + aranges[owner];
655     }
656     PetscCall(ISLocalToGlobalMappingRestoreIndices(l2gmap_p, &idxs));
657     PetscCall(ISCreateGeneral(comm, nPl, l2g_indices_p, PETSC_OWN_POINTER, &IS_l2g_p));
658 
659     /* local to global scatter for pressure */
660     PetscCall(VecScatterCreate(fetidpmat_ctx->vP, NULL, fetidp_global, IS_l2g_p, &fetidpmat_ctx->l2g_p));
661     PetscCall(ISDestroy(&IS_l2g_p));
662 
663     /* scatter for lagrange multipliers only */
664     PetscCall(VecCreate(comm, &v));
665     PetscCall(VecSetType(v, VECSTANDARD));
666     PetscCall(VecSetLayout(v, llay));
667     PetscCall(VecSetUp(v));
668     PetscCall(ISCreateGeneral(comm, n_local_lambda, l2g_indices, PETSC_COPY_VALUES, &ais));
669     PetscCall(VecScatterCreate(fetidpmat_ctx->lambda_local, NULL, v, ais, &fetidpmat_ctx->l2g_lambda_only));
670     PetscCall(ISDestroy(&ais));
671     PetscCall(VecDestroy(&v));
672 
673     /* shift local to global indices for multipliers */
674     for (i = 0; i < n_local_lambda; i++) {
675       PetscInt    ps;
676       PetscMPIInt owner;
677 
678       PetscCall(PetscLayoutFindOwner(llay, l2g_indices[i], &owner));
679       ps             = pranges[owner + 1] - pranges[owner];
680       l2g_indices[i] = l2g_indices[i] - lranges[owner] + aranges[owner] + ps;
681     }
682 
683     /* scatter from alldofs to pressures global fetidp vector */
684     PetscCall(PetscLayoutGetRange(alay, &rst, NULL));
685     PetscCall(ISCreateStride(comm, nPgl, rst, 1, &ais));
686     PetscCall(VecScatterCreate(pcis->vec1_global, pP, fetidp_global, ais, &fetidpmat_ctx->g2g_p));
687     PetscCall(ISDestroy(&ais));
688   }
689   PetscCall(PetscLayoutDestroy(&llay));
690   PetscCall(PetscLayoutDestroy(&play));
691   PetscCall(ISCreateGeneral(comm, n_local_lambda, l2g_indices, PETSC_OWN_POINTER, &IS_l2g_lambda));
692 
693   /* scatter from local to global multipliers */
694   PetscCall(VecScatterCreate(fetidpmat_ctx->lambda_local, NULL, fetidp_global, IS_l2g_lambda, &fetidpmat_ctx->l2g_lambda));
695   PetscCall(ISDestroy(&IS_l2g_lambda));
696   PetscCall(ISLocalToGlobalMappingDestroy(&l2gmap_p));
697   PetscCall(VecDestroy(&fetidp_global));
698 
699   /* Create some work vectors needed by fetidp */
700   PetscCall(VecDuplicate(pcis->vec1_B, &fetidpmat_ctx->temp_solution_B));
701   PetscCall(VecDuplicate(pcis->vec1_D, &fetidpmat_ctx->temp_solution_D));
702   PetscFunctionReturn(PETSC_SUCCESS);
703 }
704 
705 PetscErrorCode PCBDDCSetupFETIDPPCContext(Mat fetimat, FETIDPPC_ctx fetidppc_ctx)
706 {
707   FETIDPMat_ctx mat_ctx;
708   PC_BDDC      *pcbddc = (PC_BDDC *)fetidppc_ctx->pc->data;
709   PC_IS        *pcis   = (PC_IS *)fetidppc_ctx->pc->data;
710   PetscBool     lumped = PETSC_FALSE;
711 
712   PetscFunctionBegin;
713   PetscCall(MatShellGetContext(fetimat, &mat_ctx));
714   /* get references from objects created when setting up feti mat context */
715   PetscCall(PetscObjectReference((PetscObject)mat_ctx->lambda_local));
716   fetidppc_ctx->lambda_local = mat_ctx->lambda_local;
717   PetscCall(PetscObjectReference((PetscObject)mat_ctx->B_Ddelta));
718   fetidppc_ctx->B_Ddelta = mat_ctx->B_Ddelta;
719   if (mat_ctx->deluxe_nonred) {
720     PC            pc, mpc;
721     BDdelta_DN    ctx;
722     MatSolverType solver;
723     const char   *prefix;
724 
725     PetscCall(MatShellGetContext(mat_ctx->B_Ddelta, &ctx));
726     PetscCall(KSPSetType(ctx->kBD, KSPPREONLY));
727     PetscCall(KSPGetPC(ctx->kBD, &mpc));
728     PetscCall(KSPGetPC(pcbddc->ksp_D, &pc));
729     PetscCall(PCSetType(mpc, PCLU));
730     PetscCall(PCFactorGetMatSolverType(pc, &solver));
731     if (solver) PetscCall(PCFactorSetMatSolverType(mpc, solver));
732     PetscCall(MatGetOptionsPrefix(fetimat, &prefix));
733     PetscCall(KSPSetOptionsPrefix(ctx->kBD, prefix));
734     PetscCall(KSPAppendOptionsPrefix(ctx->kBD, "bddelta_"));
735     PetscCall(KSPSetFromOptions(ctx->kBD));
736   }
737 
738   if (mat_ctx->l2g_lambda_only) {
739     PetscCall(PetscObjectReference((PetscObject)mat_ctx->l2g_lambda_only));
740     fetidppc_ctx->l2g_lambda = mat_ctx->l2g_lambda_only;
741   } else {
742     PetscCall(PetscObjectReference((PetscObject)mat_ctx->l2g_lambda));
743     fetidppc_ctx->l2g_lambda = mat_ctx->l2g_lambda;
744   }
745   /* Dirichlet preconditioner */
746   PetscCall(PetscOptionsGetBool(NULL, ((PetscObject)fetimat)->prefix, "-pc_lumped", &lumped, NULL));
747   if (!lumped) {
748     IS        iV;
749     PetscBool discrete_harmonic = PETSC_FALSE;
750 
751     PetscCall(PetscObjectQuery((PetscObject)fetidppc_ctx->pc, "__KSPFETIDP_iV", (PetscObject *)&iV));
752     if (iV) PetscCall(PetscOptionsGetBool(NULL, ((PetscObject)fetimat)->prefix, "-pc_discrete_harmonic", &discrete_harmonic, NULL));
753     if (discrete_harmonic) {
754       KSP             sksp;
755       PC              pc;
756       PCBDDCSubSchurs sub_schurs = pcbddc->sub_schurs;
757       Mat             A_II, A_IB, A_BI;
758       IS              iP = NULL;
759       PetscBool       isshell, reuse = PETSC_FALSE;
760       KSPType         ksptype;
761       const char     *prefix;
762 
763       /*
764         We constructs a Schur complement for
765 
766         | A_II A_ID |
767         | A_DI A_DD |
768 
769         instead of
770 
771         | A_II  B^t_II A_ID |
772         | B_II -C_II   B_ID |
773         | A_DI  B^t_ID A_DD |
774 
775       */
776       if (sub_schurs && sub_schurs->reuse_solver) {
777         PetscCall(PetscObjectQuery((PetscObject)sub_schurs->A, "__KSPFETIDP_iP", (PetscObject *)&iP));
778         if (iP) reuse = PETSC_TRUE;
779       }
780       if (!reuse) {
781         IS       aB;
782         PetscInt nb;
783         PetscCall(ISGetLocalSize(pcis->is_B_local, &nb));
784         PetscCall(ISCreateStride(PetscObjectComm((PetscObject)pcis->A_II), nb, 0, 1, &aB));
785         PetscCall(MatCreateSubMatrix(pcis->A_II, iV, iV, MAT_INITIAL_MATRIX, &A_II));
786         PetscCall(MatCreateSubMatrix(pcis->A_IB, iV, aB, MAT_INITIAL_MATRIX, &A_IB));
787         PetscCall(MatCreateSubMatrix(pcis->A_BI, aB, iV, MAT_INITIAL_MATRIX, &A_BI));
788         PetscCall(ISDestroy(&aB));
789       } else {
790         PetscCall(MatCreateSubMatrix(sub_schurs->A, pcis->is_I_local, pcis->is_B_local, MAT_INITIAL_MATRIX, &A_IB));
791         PetscCall(MatCreateSubMatrix(sub_schurs->A, pcis->is_B_local, pcis->is_I_local, MAT_INITIAL_MATRIX, &A_BI));
792         PetscCall(PetscObjectReference((PetscObject)pcis->A_II));
793         A_II = pcis->A_II;
794       }
795       PetscCall(MatCreateSchurComplement(A_II, A_II, A_IB, A_BI, pcis->A_BB, &fetidppc_ctx->S_j));
796 
797       /* propagate settings of solver */
798       PetscCall(MatSchurComplementGetKSP(fetidppc_ctx->S_j, &sksp));
799       PetscCall(KSPGetType(pcis->ksp_D, &ksptype));
800       PetscCall(KSPSetType(sksp, ksptype));
801       PetscCall(KSPGetPC(pcis->ksp_D, &pc));
802       PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCSHELL, &isshell));
803       if (!isshell) {
804         MatSolverType solver;
805         PCType        pctype;
806 
807         PetscCall(PCGetType(pc, &pctype));
808         PetscCall(PCFactorGetMatSolverType(pc, &solver));
809         PetscCall(KSPGetPC(sksp, &pc));
810         PetscCall(PCSetType(pc, pctype));
811         if (solver) PetscCall(PCFactorSetMatSolverType(pc, solver));
812       } else {
813         PetscCall(KSPGetPC(sksp, &pc));
814         PetscCall(PCSetType(pc, PCLU));
815       }
816       PetscCall(MatDestroy(&A_II));
817       PetscCall(MatDestroy(&A_IB));
818       PetscCall(MatDestroy(&A_BI));
819       PetscCall(MatGetOptionsPrefix(fetimat, &prefix));
820       PetscCall(KSPSetOptionsPrefix(sksp, prefix));
821       PetscCall(KSPAppendOptionsPrefix(sksp, "harmonic_"));
822       PetscCall(KSPSetFromOptions(sksp));
823       if (reuse) {
824         PetscCall(KSPSetPC(sksp, sub_schurs->reuse_solver->interior_solver));
825         PetscCall(PetscObjectIncrementTabLevel((PetscObject)sub_schurs->reuse_solver->interior_solver, (PetscObject)sksp, 0));
826       }
827     } else { /* default Dirichlet preconditioner is pde-harmonic */
828       PetscCall(MatCreateSchurComplement(pcis->A_II, pcis->A_II, pcis->A_IB, pcis->A_BI, pcis->A_BB, &fetidppc_ctx->S_j));
829       PetscCall(MatSchurComplementSetKSP(fetidppc_ctx->S_j, pcis->ksp_D));
830     }
831   } else {
832     PetscCall(PetscObjectReference((PetscObject)pcis->A_BB));
833     fetidppc_ctx->S_j = pcis->A_BB;
834   }
835   /* saddle-point */
836   if (mat_ctx->xPg) {
837     PetscCall(PetscObjectReference((PetscObject)mat_ctx->xPg));
838     fetidppc_ctx->xPg = mat_ctx->xPg;
839     PetscCall(PetscObjectReference((PetscObject)mat_ctx->yPg));
840     fetidppc_ctx->yPg = mat_ctx->yPg;
841   }
842   PetscFunctionReturn(PETSC_SUCCESS);
843 }
844 
845 static PetscErrorCode FETIDPMatMult_Kernel(Mat fetimat, Vec x, Vec y, PetscBool trans)
846 {
847   FETIDPMat_ctx mat_ctx;
848   PC_BDDC      *pcbddc;
849   PC_IS        *pcis;
850 
851   PetscFunctionBegin;
852   PetscCall(MatShellGetContext(fetimat, &mat_ctx));
853   pcis   = (PC_IS *)mat_ctx->pc->data;
854   pcbddc = (PC_BDDC *)mat_ctx->pc->data;
855   /* Application of B_delta^T */
856   PetscCall(VecSet(pcis->vec1_B, 0.));
857   PetscCall(VecScatterBegin(mat_ctx->l2g_lambda, x, mat_ctx->lambda_local, INSERT_VALUES, SCATTER_REVERSE));
858   PetscCall(VecScatterEnd(mat_ctx->l2g_lambda, x, mat_ctx->lambda_local, INSERT_VALUES, SCATTER_REVERSE));
859   PetscCall(MatMultTranspose(mat_ctx->B_delta, mat_ctx->lambda_local, pcis->vec1_B));
860 
861   /* Add contribution from saddle point */
862   if (mat_ctx->l2g_p) {
863     PetscCall(VecScatterBegin(mat_ctx->l2g_p, x, mat_ctx->vP, INSERT_VALUES, SCATTER_REVERSE));
864     PetscCall(VecScatterEnd(mat_ctx->l2g_p, x, mat_ctx->vP, INSERT_VALUES, SCATTER_REVERSE));
865     if (pcbddc->switch_static) {
866       if (trans) {
867         PetscCall(MatMultTranspose(mat_ctx->B_BI, mat_ctx->vP, pcis->vec1_D));
868       } else {
869         PetscCall(MatMult(mat_ctx->Bt_BI, mat_ctx->vP, pcis->vec1_D));
870       }
871     }
872     if (trans) {
873       PetscCall(MatMultTransposeAdd(mat_ctx->B_BB, mat_ctx->vP, pcis->vec1_B, pcis->vec1_B));
874     } else {
875       PetscCall(MatMultAdd(mat_ctx->Bt_BB, mat_ctx->vP, pcis->vec1_B, pcis->vec1_B));
876     }
877   } else {
878     if (pcbddc->switch_static) PetscCall(VecSet(pcis->vec1_D, 0.0));
879   }
880   /* Application of \widetilde{S}^-1 */
881   PetscCall(PetscArrayzero(pcbddc->benign_p0, pcbddc->benign_n));
882   PetscCall(PCBDDCApplyInterfacePreconditioner(mat_ctx->pc, trans));
883   PetscCall(PetscArrayzero(pcbddc->benign_p0, pcbddc->benign_n));
884   PetscCall(VecSet(y, 0.0));
885   /* Application of B_delta */
886   PetscCall(MatMult(mat_ctx->B_delta, pcis->vec1_B, mat_ctx->lambda_local));
887   /* Contribution from boundary pressures */
888   if (mat_ctx->C) {
889     const PetscScalar *lx;
890     PetscScalar       *ly;
891 
892     /* pressure ordered first in the local part of x and y */
893     PetscCall(VecGetArrayRead(x, &lx));
894     PetscCall(VecGetArray(y, &ly));
895     PetscCall(VecPlaceArray(mat_ctx->xPg, lx));
896     PetscCall(VecPlaceArray(mat_ctx->yPg, ly));
897     if (trans) {
898       PetscCall(MatMultTranspose(mat_ctx->C, mat_ctx->xPg, mat_ctx->yPg));
899     } else {
900       PetscCall(MatMult(mat_ctx->C, mat_ctx->xPg, mat_ctx->yPg));
901     }
902     PetscCall(VecResetArray(mat_ctx->xPg));
903     PetscCall(VecResetArray(mat_ctx->yPg));
904     PetscCall(VecRestoreArrayRead(x, &lx));
905     PetscCall(VecRestoreArray(y, &ly));
906   }
907   /* Add contribution from saddle point */
908   if (mat_ctx->l2g_p) {
909     if (trans) {
910       PetscCall(MatMultTranspose(mat_ctx->Bt_BB, pcis->vec1_B, mat_ctx->vP));
911     } else {
912       PetscCall(MatMult(mat_ctx->B_BB, pcis->vec1_B, mat_ctx->vP));
913     }
914     if (pcbddc->switch_static) {
915       if (trans) {
916         PetscCall(MatMultTransposeAdd(mat_ctx->Bt_BI, pcis->vec1_D, mat_ctx->vP, mat_ctx->vP));
917       } else {
918         PetscCall(MatMultAdd(mat_ctx->B_BI, pcis->vec1_D, mat_ctx->vP, mat_ctx->vP));
919       }
920     }
921     PetscCall(VecScatterBegin(mat_ctx->l2g_p, mat_ctx->vP, y, ADD_VALUES, SCATTER_FORWARD));
922     PetscCall(VecScatterEnd(mat_ctx->l2g_p, mat_ctx->vP, y, ADD_VALUES, SCATTER_FORWARD));
923   }
924   PetscCall(VecScatterBegin(mat_ctx->l2g_lambda, mat_ctx->lambda_local, y, ADD_VALUES, SCATTER_FORWARD));
925   PetscCall(VecScatterEnd(mat_ctx->l2g_lambda, mat_ctx->lambda_local, y, ADD_VALUES, SCATTER_FORWARD));
926   PetscFunctionReturn(PETSC_SUCCESS);
927 }
928 
929 PetscErrorCode FETIDPMatMult(Mat fetimat, Vec x, Vec y)
930 {
931   PetscFunctionBegin;
932   PetscCall(FETIDPMatMult_Kernel(fetimat, x, y, PETSC_FALSE));
933   PetscFunctionReturn(PETSC_SUCCESS);
934 }
935 
936 PetscErrorCode FETIDPMatMultTranspose(Mat fetimat, Vec x, Vec y)
937 {
938   PetscFunctionBegin;
939   PetscCall(FETIDPMatMult_Kernel(fetimat, x, y, PETSC_TRUE));
940   PetscFunctionReturn(PETSC_SUCCESS);
941 }
942 
943 static PetscErrorCode FETIDPPCApply_Kernel(PC fetipc, Vec x, Vec y, PetscBool trans)
944 {
945   FETIDPPC_ctx pc_ctx;
946   PC_IS       *pcis;
947 
948   PetscFunctionBegin;
949   PetscCall(PCShellGetContext(fetipc, &pc_ctx));
950   pcis = (PC_IS *)pc_ctx->pc->data;
951   /* Application of B_Ddelta^T */
952   PetscCall(VecScatterBegin(pc_ctx->l2g_lambda, x, pc_ctx->lambda_local, INSERT_VALUES, SCATTER_REVERSE));
953   PetscCall(VecScatterEnd(pc_ctx->l2g_lambda, x, pc_ctx->lambda_local, INSERT_VALUES, SCATTER_REVERSE));
954   PetscCall(VecSet(pcis->vec2_B, 0.0));
955   PetscCall(MatMultTranspose(pc_ctx->B_Ddelta, pc_ctx->lambda_local, pcis->vec2_B));
956   /* Application of local Schur complement */
957   if (trans) {
958     PetscCall(MatMultTranspose(pc_ctx->S_j, pcis->vec2_B, pcis->vec1_B));
959   } else {
960     PetscCall(MatMult(pc_ctx->S_j, pcis->vec2_B, pcis->vec1_B));
961   }
962   /* Application of B_Ddelta */
963   PetscCall(MatMult(pc_ctx->B_Ddelta, pcis->vec1_B, pc_ctx->lambda_local));
964   PetscCall(VecSet(y, 0.0));
965   PetscCall(VecScatterBegin(pc_ctx->l2g_lambda, pc_ctx->lambda_local, y, ADD_VALUES, SCATTER_FORWARD));
966   PetscCall(VecScatterEnd(pc_ctx->l2g_lambda, pc_ctx->lambda_local, y, ADD_VALUES, SCATTER_FORWARD));
967   PetscFunctionReturn(PETSC_SUCCESS);
968 }
969 
970 PetscErrorCode FETIDPPCApply(PC pc, Vec x, Vec y)
971 {
972   PetscFunctionBegin;
973   PetscCall(FETIDPPCApply_Kernel(pc, x, y, PETSC_FALSE));
974   PetscFunctionReturn(PETSC_SUCCESS);
975 }
976 
977 PetscErrorCode FETIDPPCApplyTranspose(PC pc, Vec x, Vec y)
978 {
979   PetscFunctionBegin;
980   PetscCall(FETIDPPCApply_Kernel(pc, x, y, PETSC_TRUE));
981   PetscFunctionReturn(PETSC_SUCCESS);
982 }
983 
984 PetscErrorCode FETIDPPCView(PC pc, PetscViewer viewer)
985 {
986   FETIDPPC_ctx pc_ctx;
987   PetscBool    iascii;
988   PetscViewer  sviewer;
989 
990   PetscFunctionBegin;
991   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
992   if (iascii) {
993     PetscMPIInt rank;
994     PetscBool   isschur, isshell;
995 
996     PetscCall(PCShellGetContext(pc, &pc_ctx));
997     PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank));
998     PetscCall(PetscObjectTypeCompare((PetscObject)pc_ctx->S_j, MATSCHURCOMPLEMENT, &isschur));
999     if (isschur) {
1000       PetscCall(PetscViewerASCIIPrintf(viewer, "  Dirichlet preconditioner (just from rank 0)\n"));
1001     } else {
1002       PetscCall(PetscViewerASCIIPrintf(viewer, "  Lumped preconditioner (just from rank 0)\n"));
1003     }
1004     PetscCall(PetscViewerGetSubViewer(viewer, PetscObjectComm((PetscObject)pc_ctx->S_j), &sviewer));
1005     if (rank == 0) {
1006       PetscCall(PetscViewerPushFormat(sviewer, PETSC_VIEWER_ASCII_INFO));
1007       PetscCall(PetscViewerASCIIPushTab(sviewer));
1008       PetscCall(MatView(pc_ctx->S_j, sviewer));
1009       PetscCall(PetscViewerASCIIPopTab(sviewer));
1010       PetscCall(PetscViewerPopFormat(sviewer));
1011     }
1012     PetscCall(PetscViewerRestoreSubViewer(viewer, PetscObjectComm((PetscObject)pc_ctx->S_j), &sviewer));
1013     PetscCall(PetscObjectTypeCompare((PetscObject)pc_ctx->B_Ddelta, MATSHELL, &isshell));
1014     if (isshell) {
1015       BDdelta_DN ctx;
1016       PetscCall(PetscViewerASCIIPrintf(viewer, "  FETI-DP BDdelta: DB^t * (B D^-1 B^t)^-1 for deluxe scaling (just from rank 0)\n"));
1017       PetscCall(MatShellGetContext(pc_ctx->B_Ddelta, &ctx));
1018       PetscCall(PetscViewerGetSubViewer(viewer, PetscObjectComm((PetscObject)pc_ctx->S_j), &sviewer));
1019       if (rank == 0) {
1020         PetscInt tl;
1021 
1022         PetscCall(PetscViewerASCIIGetTab(sviewer, &tl));
1023         PetscCall(PetscObjectSetTabLevel((PetscObject)ctx->kBD, tl));
1024         PetscCall(KSPView(ctx->kBD, sviewer));
1025         PetscCall(PetscViewerPushFormat(sviewer, PETSC_VIEWER_ASCII_INFO));
1026         PetscCall(MatView(ctx->BD, sviewer));
1027         PetscCall(PetscViewerPopFormat(sviewer));
1028       }
1029       PetscCall(PetscViewerRestoreSubViewer(viewer, PetscObjectComm((PetscObject)pc_ctx->S_j), &sviewer));
1030     }
1031     PetscCall(PetscViewerFlush(viewer));
1032   }
1033   PetscFunctionReturn(PETSC_SUCCESS);
1034 }
1035