xref: /petsc/src/ksp/pc/impls/bddc/bddcfetidp.c (revision 9de2952e64edfa51c4ec5537730b8a834bdd3686)
15e5bbd0aSStefano Zampini #include <petsc/private/pcbddcimpl.h>
25e5bbd0aSStefano Zampini #include <petsc/private/pcbddcprivateimpl.h>
325d06dbeSstefano_zampini #include <petscblaslapack.h>
425d06dbeSstefano_zampini 
5d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMult_BDdelta_deluxe_nonred(Mat A, Vec x, Vec y)
6d71ae5a4SJacob Faibussowitsch {
725d06dbeSstefano_zampini   BDdelta_DN ctx;
825d06dbeSstefano_zampini 
925d06dbeSstefano_zampini   PetscFunctionBegin;
109566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(A, &ctx));
119566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(ctx->BD, x, ctx->work));
129566063dSJacob Faibussowitsch   PetscCall(KSPSolveTranspose(ctx->kBD, ctx->work, y));
13c0decd05SBarry Smith   /* No PC so cannot propagate up failure in KSPSolveTranspose() */
143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1525d06dbeSstefano_zampini }
1625d06dbeSstefano_zampini 
17d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMultTranspose_BDdelta_deluxe_nonred(Mat A, Vec x, Vec y)
18d71ae5a4SJacob Faibussowitsch {
1925d06dbeSstefano_zampini   BDdelta_DN ctx;
2025d06dbeSstefano_zampini 
2125d06dbeSstefano_zampini   PetscFunctionBegin;
229566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(A, &ctx));
239566063dSJacob Faibussowitsch   PetscCall(KSPSolve(ctx->kBD, x, ctx->work));
24c0decd05SBarry Smith   /* No PC so cannot propagate up failure in KSPSolve() */
259566063dSJacob Faibussowitsch   PetscCall(MatMult(ctx->BD, ctx->work, y));
263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2725d06dbeSstefano_zampini }
2825d06dbeSstefano_zampini 
29d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatDestroy_BDdelta_deluxe_nonred(Mat A)
30d71ae5a4SJacob Faibussowitsch {
3125d06dbeSstefano_zampini   BDdelta_DN ctx;
3225d06dbeSstefano_zampini 
3325d06dbeSstefano_zampini   PetscFunctionBegin;
349566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(A, &ctx));
359566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&ctx->BD));
369566063dSJacob Faibussowitsch   PetscCall(KSPDestroy(&ctx->kBD));
379566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ctx->work));
389566063dSJacob Faibussowitsch   PetscCall(PetscFree(ctx));
393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4025d06dbeSstefano_zampini }
4125d06dbeSstefano_zampini 
42d71ae5a4SJacob Faibussowitsch PetscErrorCode PCBDDCCreateFETIDPMatContext(PC pc, FETIDPMat_ctx *fetidpmat_ctx)
43d71ae5a4SJacob Faibussowitsch {
44674ae819SStefano Zampini   FETIDPMat_ctx newctx;
45674ae819SStefano Zampini 
46674ae819SStefano Zampini   PetscFunctionBegin;
479566063dSJacob Faibussowitsch   PetscCall(PetscNew(&newctx));
48674ae819SStefano Zampini   /* increase the reference count for BDDC preconditioner */
499566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)pc));
50674ae819SStefano Zampini   newctx->pc     = pc;
51674ae819SStefano Zampini   *fetidpmat_ctx = newctx;
523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
53674ae819SStefano Zampini }
54674ae819SStefano Zampini 
55d71ae5a4SJacob Faibussowitsch PetscErrorCode PCBDDCCreateFETIDPPCContext(PC pc, FETIDPPC_ctx *fetidppc_ctx)
56d71ae5a4SJacob Faibussowitsch {
57674ae819SStefano Zampini   FETIDPPC_ctx newctx;
58674ae819SStefano Zampini 
59674ae819SStefano Zampini   PetscFunctionBegin;
609566063dSJacob Faibussowitsch   PetscCall(PetscNew(&newctx));
61674ae819SStefano Zampini   /* increase the reference count for BDDC preconditioner */
629566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)pc));
63674ae819SStefano Zampini   newctx->pc    = pc;
64674ae819SStefano Zampini   *fetidppc_ctx = newctx;
653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
66674ae819SStefano Zampini }
67674ae819SStefano Zampini 
68d71ae5a4SJacob Faibussowitsch PetscErrorCode PCBDDCDestroyFETIDPMat(Mat A)
69d71ae5a4SJacob Faibussowitsch {
70674ae819SStefano Zampini   FETIDPMat_ctx mat_ctx;
71674ae819SStefano Zampini 
72674ae819SStefano Zampini   PetscFunctionBegin;
739566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(A, &mat_ctx));
749566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&mat_ctx->lambda_local));
759566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&mat_ctx->temp_solution_D));
769566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&mat_ctx->temp_solution_B));
779566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&mat_ctx->B_delta));
789566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&mat_ctx->B_Ddelta));
799566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&mat_ctx->B_BB));
809566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&mat_ctx->B_BI));
819566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&mat_ctx->Bt_BB));
829566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&mat_ctx->Bt_BI));
839566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&mat_ctx->C));
849566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&mat_ctx->rhs_flip));
859566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&mat_ctx->vP));
869566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&mat_ctx->xPg));
879566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&mat_ctx->yPg));
889566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&mat_ctx->l2g_lambda));
899566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&mat_ctx->l2g_lambda_only));
909566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&mat_ctx->l2g_p));
919566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&mat_ctx->g2g_p));
929566063dSJacob Faibussowitsch   PetscCall(PCDestroy(&mat_ctx->pc)); /* decrease PCBDDC reference count */
939566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&mat_ctx->pressure));
949566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&mat_ctx->lagrange));
959566063dSJacob Faibussowitsch   PetscCall(PetscFree(mat_ctx));
963ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
97674ae819SStefano Zampini }
98674ae819SStefano Zampini 
99d71ae5a4SJacob Faibussowitsch PetscErrorCode PCBDDCDestroyFETIDPPC(PC pc)
100d71ae5a4SJacob Faibussowitsch {
101674ae819SStefano Zampini   FETIDPPC_ctx pc_ctx;
102674ae819SStefano Zampini 
103674ae819SStefano Zampini   PetscFunctionBegin;
1049566063dSJacob Faibussowitsch   PetscCall(PCShellGetContext(pc, &pc_ctx));
1059566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pc_ctx->lambda_local));
1069566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pc_ctx->B_Ddelta));
1079566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&pc_ctx->l2g_lambda));
1089566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pc_ctx->S_j));
1099566063dSJacob Faibussowitsch   PetscCall(PCDestroy(&pc_ctx->pc)); /* decrease PCBDDC reference count */
1109566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pc_ctx->xPg));
1119566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pc_ctx->yPg));
1129566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc_ctx));
1133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
114674ae819SStefano Zampini }
115674ae819SStefano Zampini 
116d71ae5a4SJacob Faibussowitsch PetscErrorCode PCBDDCSetupFETIDPMatContext(FETIDPMat_ctx fetidpmat_ctx)
117d71ae5a4SJacob Faibussowitsch {
118674ae819SStefano Zampini   PC_IS          *pcis      = (PC_IS *)fetidpmat_ctx->pc->data;
119674ae819SStefano Zampini   PC_BDDC        *pcbddc    = (PC_BDDC *)fetidpmat_ctx->pc->data;
120674ae819SStefano Zampini   PCBDDCGraph     mat_graph = pcbddc->mat_graph;
121674ae819SStefano Zampini   Mat_IS         *matis     = (Mat_IS *)fetidpmat_ctx->pc->pmat->data;
122674ae819SStefano Zampini   MPI_Comm        comm;
12325d06dbeSstefano_zampini   Mat             ScalingMat, BD1, BD2;
124457d4a33Sstefano_zampini   Vec             fetidp_global;
125674ae819SStefano Zampini   IS              IS_l2g_lambda;
126a1c0d0daSstefano_zampini   IS              subset, subset_mult, subset_n, isvert;
127674ae819SStefano Zampini   PetscBool       skip_node, fully_redundant;
128674ae819SStefano Zampini   PetscInt        i, j, k, s, n_boundary_dofs, n_global_lambda, n_vertices, partial_sum;
129dc456d91SStefano Zampini   PetscInt        cum, n_local_lambda, n_lambda_for_dof, dual_size, n_neg_values, n_pos_values;
13076ec1555SBarry Smith   PetscMPIInt     rank, size, buf_size, neigh;
131674ae819SStefano Zampini   PetscScalar     scalar_value;
132a1c0d0daSstefano_zampini   const PetscInt *vertex_indices;
133dc456d91SStefano Zampini   PetscInt       *dual_dofs_boundary_indices, *aux_local_numbering_1;
134dc456d91SStefano Zampini   const PetscInt *aux_global_numbering;
135674ae819SStefano Zampini   PetscInt       *aux_sums, *cols_B_delta, *l2g_indices;
136674ae819SStefano Zampini   PetscScalar    *array, *scaling_factors, *vals_B_delta;
13725d06dbeSstefano_zampini   PetscScalar   **all_factors;
138674ae819SStefano Zampini   PetscInt       *aux_local_numbering_2;
139*9de2952eSStefano Zampini   PetscInt       *count, **neighbours_set;
140457d4a33Sstefano_zampini   PetscLayout     llay;
141a1c0d0daSstefano_zampini 
142457d4a33Sstefano_zampini   /* saddlepoint */
143457d4a33Sstefano_zampini   ISLocalToGlobalMapping l2gmap_p;
144457d4a33Sstefano_zampini   PetscLayout            play;
145457d4a33Sstefano_zampini   IS                     gP, pP;
146457d4a33Sstefano_zampini   PetscInt               nPl, nPg, nPgl;
147674ae819SStefano Zampini 
148674ae819SStefano Zampini   PetscFunctionBegin;
149f4f49eeaSPierre Jolivet   PetscCall(PetscObjectGetComm((PetscObject)fetidpmat_ctx->pc, &comm));
1509566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
1519566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
152674ae819SStefano Zampini 
153457d4a33Sstefano_zampini   /* saddlepoint */
154457d4a33Sstefano_zampini   nPl      = 0;
155457d4a33Sstefano_zampini   nPg      = 0;
156457d4a33Sstefano_zampini   nPgl     = 0;
157457d4a33Sstefano_zampini   gP       = NULL;
158457d4a33Sstefano_zampini   pP       = NULL;
159457d4a33Sstefano_zampini   l2gmap_p = NULL;
160457d4a33Sstefano_zampini   play     = NULL;
1619566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)fetidpmat_ctx->pc, "__KSPFETIDP_pP", (PetscObject *)&pP));
162022d8d2bSstefano_zampini   if (pP) { /* saddle point */
163457d4a33Sstefano_zampini     /* subdomain pressures in global numbering */
1649566063dSJacob Faibussowitsch     PetscCall(PetscObjectQuery((PetscObject)fetidpmat_ctx->pc, "__KSPFETIDP_gP", (PetscObject *)&gP));
16528b400f6SJacob Faibussowitsch     PetscCheck(gP, PETSC_COMM_SELF, PETSC_ERR_PLIB, "gP not present");
1669566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(gP, &nPl));
1679566063dSJacob Faibussowitsch     PetscCall(VecCreate(PETSC_COMM_SELF, &fetidpmat_ctx->vP));
1689566063dSJacob Faibussowitsch     PetscCall(VecSetSizes(fetidpmat_ctx->vP, nPl, nPl));
1699566063dSJacob Faibussowitsch     PetscCall(VecSetType(fetidpmat_ctx->vP, VECSTANDARD));
1709566063dSJacob Faibussowitsch     PetscCall(VecSetUp(fetidpmat_ctx->vP));
171457d4a33Sstefano_zampini 
172e1214c54Sstefano_zampini     /* pressure matrix */
1739566063dSJacob Faibussowitsch     PetscCall(PetscObjectQuery((PetscObject)fetidpmat_ctx->pc, "__KSPFETIDP_C", (PetscObject *)&fetidpmat_ctx->C));
174e1214c54Sstefano_zampini     if (!fetidpmat_ctx->C) { /* null pressure block, compute layout and global numbering for pressures */
175457d4a33Sstefano_zampini       IS Pg;
176457d4a33Sstefano_zampini 
1779566063dSJacob Faibussowitsch       PetscCall(ISRenumber(gP, NULL, &nPg, &Pg));
1789566063dSJacob Faibussowitsch       PetscCall(ISLocalToGlobalMappingCreateIS(Pg, &l2gmap_p));
1799566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&Pg));
1809566063dSJacob Faibussowitsch       PetscCall(PetscLayoutCreate(comm, &play));
1819566063dSJacob Faibussowitsch       PetscCall(PetscLayoutSetBlockSize(play, 1));
1829566063dSJacob Faibussowitsch       PetscCall(PetscLayoutSetSize(play, nPg));
1839566063dSJacob Faibussowitsch       PetscCall(ISGetLocalSize(pP, &nPgl));
1849566063dSJacob Faibussowitsch       PetscCall(PetscLayoutSetLocalSize(play, nPgl));
1859566063dSJacob Faibussowitsch       PetscCall(PetscLayoutSetUp(play));
186457d4a33Sstefano_zampini     } else {
1879566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)fetidpmat_ctx->C));
1889566063dSJacob Faibussowitsch       PetscCall(MatISGetLocalToGlobalMapping(fetidpmat_ctx->C, &l2gmap_p, NULL));
1899566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)l2gmap_p));
1909566063dSJacob Faibussowitsch       PetscCall(MatGetSize(fetidpmat_ctx->C, &nPg, NULL));
1919566063dSJacob Faibussowitsch       PetscCall(MatGetLocalSize(fetidpmat_ctx->C, NULL, &nPgl));
1929566063dSJacob Faibussowitsch       PetscCall(MatGetLayouts(fetidpmat_ctx->C, NULL, &llay));
1939566063dSJacob Faibussowitsch       PetscCall(PetscLayoutReference(llay, &play));
194457d4a33Sstefano_zampini     }
1959566063dSJacob Faibussowitsch     PetscCall(VecCreateMPIWithArray(comm, 1, nPgl, nPg, NULL, &fetidpmat_ctx->xPg));
1969566063dSJacob Faibussowitsch     PetscCall(VecCreateMPIWithArray(comm, 1, nPgl, nPg, NULL, &fetidpmat_ctx->yPg));
197457d4a33Sstefano_zampini 
198e1214c54Sstefano_zampini     /* import matrices for pressures coupling */
1999566063dSJacob Faibussowitsch     PetscCall(PetscObjectQuery((PetscObject)fetidpmat_ctx->pc, "__KSPFETIDP_B_BI", (PetscObject *)&fetidpmat_ctx->B_BI));
20028b400f6SJacob Faibussowitsch     PetscCheck(fetidpmat_ctx->B_BI, PETSC_COMM_SELF, PETSC_ERR_PLIB, "B_BI not present");
2019566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)fetidpmat_ctx->B_BI));
202457d4a33Sstefano_zampini 
2039566063dSJacob Faibussowitsch     PetscCall(PetscObjectQuery((PetscObject)fetidpmat_ctx->pc, "__KSPFETIDP_B_BB", (PetscObject *)&fetidpmat_ctx->B_BB));
20428b400f6SJacob Faibussowitsch     PetscCheck(fetidpmat_ctx->B_BB, PETSC_COMM_SELF, PETSC_ERR_PLIB, "B_BB not present");
2059566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)fetidpmat_ctx->B_BB));
206457d4a33Sstefano_zampini 
2079566063dSJacob Faibussowitsch     PetscCall(PetscObjectQuery((PetscObject)fetidpmat_ctx->pc, "__KSPFETIDP_Bt_BI", (PetscObject *)&fetidpmat_ctx->Bt_BI));
20828b400f6SJacob Faibussowitsch     PetscCheck(fetidpmat_ctx->Bt_BI, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Bt_BI not present");
2099566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)fetidpmat_ctx->Bt_BI));
210457d4a33Sstefano_zampini 
2119566063dSJacob Faibussowitsch     PetscCall(PetscObjectQuery((PetscObject)fetidpmat_ctx->pc, "__KSPFETIDP_Bt_BB", (PetscObject *)&fetidpmat_ctx->Bt_BB));
21228b400f6SJacob Faibussowitsch     PetscCheck(fetidpmat_ctx->Bt_BB, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Bt_BB not present");
2139566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)fetidpmat_ctx->Bt_BB));
2146cc1294bSstefano_zampini 
2159566063dSJacob Faibussowitsch     PetscCall(PetscObjectQuery((PetscObject)fetidpmat_ctx->pc, "__KSPFETIDP_flip", (PetscObject *)&fetidpmat_ctx->rhs_flip));
2161baa6e33SBarry Smith     if (fetidpmat_ctx->rhs_flip) PetscCall(PetscObjectReference((PetscObject)fetidpmat_ctx->rhs_flip));
217457d4a33Sstefano_zampini   }
218457d4a33Sstefano_zampini 
219674ae819SStefano Zampini   /* Default type of lagrange multipliers is non-redundant */
220329cd39dSStefano Zampini   fully_redundant = fetidpmat_ctx->fully_redundant;
221674ae819SStefano Zampini 
222674ae819SStefano Zampini   /* Evaluate local and global number of lagrange multipliers */
2239566063dSJacob Faibussowitsch   PetscCall(VecSet(pcis->vec1_N, 0.0));
224674ae819SStefano Zampini   n_local_lambda  = 0;
225674ae819SStefano Zampini   partial_sum     = 0;
226674ae819SStefano Zampini   n_boundary_dofs = 0;
227a1c0d0daSstefano_zampini 
228674ae819SStefano Zampini   /* Get Vertices used to define the BDDC */
229*9de2952eSStefano Zampini   PetscCall(PCBDDCGraphGetCandidatesIS(mat_graph, NULL, NULL, NULL, NULL, &isvert));
2309566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(isvert, &n_vertices));
2319566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(isvert, &vertex_indices));
232a1c0d0daSstefano_zampini 
233674ae819SStefano Zampini   dual_size = pcis->n_B - n_vertices;
2349566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(dual_size, &dual_dofs_boundary_indices));
2359566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(dual_size, &aux_local_numbering_1));
2369566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(dual_size, &aux_local_numbering_2));
237674ae819SStefano Zampini 
238*9de2952eSStefano Zampini   /* the code below does not support multiple subdomains per process
239*9de2952eSStefano Zampini      error out in this case
240*9de2952eSStefano Zampini      TODO: I guess I can use PetscSFGetMultiSF and the code will be easier and more general */
241*9de2952eSStefano Zampini   PetscCall(PetscMalloc2(pcis->n, &count, pcis->n, &neighbours_set));
242*9de2952eSStefano Zampini   for (i = 0, j = 0; i < pcis->n; i++) j += mat_graph->nodes[i].count;
243*9de2952eSStefano Zampini   if (pcis->n) PetscCall(PetscMalloc1(j, &neighbours_set[0]));
244674ae819SStefano Zampini   for (i = 0; i < pcis->n; i++) {
245*9de2952eSStefano Zampini     PCBDDCGraphNode *node = &mat_graph->nodes[i];
246*9de2952eSStefano Zampini 
247*9de2952eSStefano Zampini     count[i] = 0;
248*9de2952eSStefano Zampini     for (j = 0; j < node->count; j++) {
249*9de2952eSStefano Zampini       if (node->neighbours_set[j] == rank) continue;
250*9de2952eSStefano Zampini       neighbours_set[i][count[i]++] = node->neighbours_set[j];
251*9de2952eSStefano Zampini     }
252*9de2952eSStefano Zampini     PetscCheck(count[i] == node->count - 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Multiple subdomains per process not supported");
253*9de2952eSStefano Zampini     s = count[i];
254*9de2952eSStefano Zampini     PetscCall(PetscSortRemoveDupsInt(count + i, neighbours_set[i]));
255*9de2952eSStefano Zampini     PetscCheck(s == count[i], PETSC_COMM_SELF, PETSC_ERR_SUP, "Multiple subdomains per process not supported");
256*9de2952eSStefano Zampini     if (i != pcis->n - 1) neighbours_set[i + 1] = neighbours_set[i] + count[i];
257*9de2952eSStefano Zampini   }
258*9de2952eSStefano Zampini 
259*9de2952eSStefano Zampini   PetscCall(VecGetArray(pcis->vec1_N, &array));
260*9de2952eSStefano Zampini   for (i = 0, s = 0; i < pcis->n; i++) {
261*9de2952eSStefano Zampini     j = count[i]; /* RECALL: count[i] does not count myself */
2622d456bbaSstefano_zampini     if (j > 0) n_boundary_dofs++;
263674ae819SStefano Zampini     skip_node = PETSC_FALSE;
264674ae819SStefano Zampini     if (s < n_vertices && vertex_indices[s] == i) { /* it works for a sorted set of vertices */
265674ae819SStefano Zampini       skip_node = PETSC_TRUE;
266674ae819SStefano Zampini       s++;
267674ae819SStefano Zampini     }
2682d456bbaSstefano_zampini     if (j < 1) skip_node = PETSC_TRUE;
269*9de2952eSStefano Zampini     if (mat_graph->nodes[i].special_dof == PCBDDCGRAPH_DIRICHLET_MARK) skip_node = PETSC_TRUE;
270674ae819SStefano Zampini     if (!skip_node) {
271674ae819SStefano Zampini       if (fully_redundant) {
272674ae819SStefano Zampini         /* fully redundant set of lagrange multipliers */
273674ae819SStefano Zampini         n_lambda_for_dof = (j * (j + 1)) / 2;
274674ae819SStefano Zampini       } else {
275674ae819SStefano Zampini         n_lambda_for_dof = j;
276674ae819SStefano Zampini       }
277674ae819SStefano Zampini       n_local_lambda += j;
278674ae819SStefano Zampini       /* needed to evaluate global number of lagrange multipliers */
279674ae819SStefano Zampini       array[i] = (1.0 * n_lambda_for_dof) / (j + 1.0); /* already scaled for the next global sum */
280674ae819SStefano Zampini       /* store some data needed */
281674ae819SStefano Zampini       dual_dofs_boundary_indices[partial_sum] = n_boundary_dofs - 1;
282674ae819SStefano Zampini       aux_local_numbering_1[partial_sum]      = i;
283674ae819SStefano Zampini       aux_local_numbering_2[partial_sum]      = n_lambda_for_dof;
284674ae819SStefano Zampini       partial_sum++;
285674ae819SStefano Zampini     }
286674ae819SStefano Zampini   }
2879566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(pcis->vec1_N, &array));
2889566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isvert, &vertex_indices));
289*9de2952eSStefano Zampini   PetscCall(PCBDDCGraphRestoreCandidatesIS(mat_graph, NULL, NULL, NULL, NULL, &isvert));
2902d456bbaSstefano_zampini   dual_size = partial_sum;
291674ae819SStefano Zampini 
292674ae819SStefano Zampini   /* compute global ordering of lagrange multipliers and associate l2g map */
2939566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(comm, partial_sum, aux_local_numbering_1, PETSC_COPY_VALUES, &subset_n));
2949566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingApplyIS(pcis->mapping, subset_n, &subset));
2959566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&subset_n));
2969566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(comm, partial_sum, aux_local_numbering_2, PETSC_OWN_POINTER, &subset_mult));
2979566063dSJacob Faibussowitsch   PetscCall(ISRenumber(subset, subset_mult, &fetidpmat_ctx->n_lambda, &subset_n));
2989566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&subset));
2993d996552SStefano Zampini 
30076bd3646SJed Brown   if (PetscDefined(USE_DEBUG)) {
3019566063dSJacob Faibussowitsch     PetscCall(VecSet(pcis->vec1_global, 0.0));
3029566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(matis->rctx, pcis->vec1_N, pcis->vec1_global, ADD_VALUES, SCATTER_REVERSE));
3039566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(matis->rctx, pcis->vec1_N, pcis->vec1_global, ADD_VALUES, SCATTER_REVERSE));
3049566063dSJacob Faibussowitsch     PetscCall(VecSum(pcis->vec1_global, &scalar_value));
3054fe826edSStefano Zampini     i = (PetscInt)PetscRealPart(scalar_value);
30663a3b9bcSJacob Faibussowitsch     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);
30776bd3646SJed Brown   }
308674ae819SStefano Zampini 
309674ae819SStefano Zampini   /* init data for scaling factors exchange */
31025d06dbeSstefano_zampini   if (!pcbddc->use_deluxe_scaling) {
31125d06dbeSstefano_zampini     PetscInt    *ptrs_buffer, neigh_position;
31225d06dbeSstefano_zampini     PetscScalar *send_buffer, *recv_buffer;
31325d06dbeSstefano_zampini     MPI_Request *send_reqs, *recv_reqs;
31425d06dbeSstefano_zampini 
315674ae819SStefano Zampini     partial_sum = 0;
3169566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(pcis->n_neigh, &ptrs_buffer));
3179566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(PetscMax(pcis->n_neigh - 1, 0), &send_reqs));
3189566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(PetscMax(pcis->n_neigh - 1, 0), &recv_reqs));
3199566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(pcis->n + 1, &all_factors));
3204b2aedd3SStefano Zampini     if (pcis->n_neigh > 0) ptrs_buffer[0] = 0;
321674ae819SStefano Zampini     for (i = 1; i < pcis->n_neigh; i++) {
322674ae819SStefano Zampini       partial_sum += pcis->n_shared[i];
323674ae819SStefano Zampini       ptrs_buffer[i] = ptrs_buffer[i - 1] + pcis->n_shared[i];
324674ae819SStefano Zampini     }
3259566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(partial_sum, &send_buffer));
3269566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(partial_sum, &recv_buffer));
3279566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(partial_sum, &all_factors[0]));
328674ae819SStefano Zampini     for (i = 0; i < pcis->n - 1; i++) {
329*9de2952eSStefano Zampini       j                  = count[i];
330674ae819SStefano Zampini       all_factors[i + 1] = all_factors[i] + j;
331674ae819SStefano Zampini     }
33225d06dbeSstefano_zampini 
333674ae819SStefano Zampini     /* scatter B scaling to N vec */
3349566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(pcis->N_to_B, pcis->D, pcis->vec1_N, INSERT_VALUES, SCATTER_REVERSE));
3359566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(pcis->N_to_B, pcis->D, pcis->vec1_N, INSERT_VALUES, SCATTER_REVERSE));
336674ae819SStefano Zampini     /* communications */
3379566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(pcis->vec1_N, (const PetscScalar **)&array));
338674ae819SStefano Zampini     for (i = 1; i < pcis->n_neigh; i++) {
339ad540459SPierre Jolivet       for (j = 0; j < pcis->n_shared[i]; j++) send_buffer[ptrs_buffer[i - 1] + j] = array[pcis->shared[i][j]];
3409566063dSJacob Faibussowitsch       PetscCall(PetscMPIIntCast(ptrs_buffer[i] - ptrs_buffer[i - 1], &buf_size));
3419566063dSJacob Faibussowitsch       PetscCall(PetscMPIIntCast(pcis->neigh[i], &neigh));
3429566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Isend(&send_buffer[ptrs_buffer[i - 1]], buf_size, MPIU_SCALAR, neigh, 0, comm, &send_reqs[i - 1]));
3439566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Irecv(&recv_buffer[ptrs_buffer[i - 1]], buf_size, MPIU_SCALAR, neigh, 0, comm, &recv_reqs[i - 1]));
344674ae819SStefano Zampini     }
3459566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(pcis->vec1_N, (const PetscScalar **)&array));
34648a46eb9SPierre Jolivet     if (pcis->n_neigh > 0) PetscCallMPI(MPI_Waitall(pcis->n_neigh - 1, recv_reqs, MPI_STATUSES_IGNORE));
347674ae819SStefano Zampini     /* put values in correct places */
348674ae819SStefano Zampini     for (i = 1; i < pcis->n_neigh; i++) {
349674ae819SStefano Zampini       for (j = 0; j < pcis->n_shared[i]; j++) {
350674ae819SStefano Zampini         k              = pcis->shared[i][j];
351674ae819SStefano Zampini         neigh_position = 0;
352*9de2952eSStefano Zampini         while (neighbours_set[k][neigh_position] != pcis->neigh[i]) { neigh_position++; }
353674ae819SStefano Zampini         all_factors[k][neigh_position] = recv_buffer[ptrs_buffer[i - 1] + j];
354674ae819SStefano Zampini       }
355674ae819SStefano Zampini     }
35648a46eb9SPierre Jolivet     if (pcis->n_neigh > 0) PetscCallMPI(MPI_Waitall(pcis->n_neigh - 1, send_reqs, MPI_STATUSES_IGNORE));
3579566063dSJacob Faibussowitsch     PetscCall(PetscFree(send_reqs));
3589566063dSJacob Faibussowitsch     PetscCall(PetscFree(recv_reqs));
3599566063dSJacob Faibussowitsch     PetscCall(PetscFree(send_buffer));
3609566063dSJacob Faibussowitsch     PetscCall(PetscFree(recv_buffer));
3619566063dSJacob Faibussowitsch     PetscCall(PetscFree(ptrs_buffer));
36225d06dbeSstefano_zampini   }
363674ae819SStefano Zampini 
364674ae819SStefano Zampini   /* Compute B and B_delta (local actions) */
3659566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(pcis->n_neigh, &aux_sums));
3669566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n_local_lambda, &l2g_indices));
3679566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n_local_lambda, &vals_B_delta));
3689566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n_local_lambda, &cols_B_delta));
36925d06dbeSstefano_zampini   if (!pcbddc->use_deluxe_scaling) {
3709566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n_local_lambda, &scaling_factors));
37125d06dbeSstefano_zampini   } else {
37225d06dbeSstefano_zampini     scaling_factors = NULL;
37325d06dbeSstefano_zampini     all_factors     = NULL;
37425d06dbeSstefano_zampini   }
3759566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(subset_n, &aux_global_numbering));
376674ae819SStefano Zampini   partial_sum = 0;
377dc456d91SStefano Zampini   cum         = 0;
378674ae819SStefano Zampini   for (i = 0; i < dual_size; i++) {
379dc456d91SStefano Zampini     n_global_lambda = aux_global_numbering[cum];
380*9de2952eSStefano Zampini     j               = count[aux_local_numbering_1[i]];
381674ae819SStefano Zampini     aux_sums[0]     = 0;
382ad540459SPierre Jolivet     for (s = 1; s < j; s++) aux_sums[s] = aux_sums[s - 1] + j - s + 1;
38325d06dbeSstefano_zampini     if (all_factors) array = all_factors[aux_local_numbering_1[i]];
384674ae819SStefano Zampini     n_neg_values = 0;
385*9de2952eSStefano Zampini     while (n_neg_values < j && neighbours_set[aux_local_numbering_1[i]][n_neg_values] < rank) { n_neg_values++; }
386674ae819SStefano Zampini     n_pos_values = j - n_neg_values;
387674ae819SStefano Zampini     if (fully_redundant) {
388674ae819SStefano Zampini       for (s = 0; s < n_neg_values; s++) {
389674ae819SStefano Zampini         l2g_indices[partial_sum + s]  = aux_sums[s] + n_neg_values - s - 1 + n_global_lambda;
390674ae819SStefano Zampini         cols_B_delta[partial_sum + s] = dual_dofs_boundary_indices[i];
391674ae819SStefano Zampini         vals_B_delta[partial_sum + s] = -1.0;
39225d06dbeSstefano_zampini         if (!pcbddc->use_deluxe_scaling) scaling_factors[partial_sum + s] = array[s];
393674ae819SStefano Zampini       }
394674ae819SStefano Zampini       for (s = 0; s < n_pos_values; s++) {
395674ae819SStefano Zampini         l2g_indices[partial_sum + s + n_neg_values]  = aux_sums[n_neg_values] + s + n_global_lambda;
396674ae819SStefano Zampini         cols_B_delta[partial_sum + s + n_neg_values] = dual_dofs_boundary_indices[i];
397674ae819SStefano Zampini         vals_B_delta[partial_sum + s + n_neg_values] = 1.0;
39825d06dbeSstefano_zampini         if (!pcbddc->use_deluxe_scaling) scaling_factors[partial_sum + s + n_neg_values] = array[s + n_neg_values];
399674ae819SStefano Zampini       }
400674ae819SStefano Zampini       partial_sum += j;
401674ae819SStefano Zampini     } else {
402674ae819SStefano Zampini       /* l2g_indices and default cols and vals of B_delta */
403674ae819SStefano Zampini       for (s = 0; s < j; s++) {
404674ae819SStefano Zampini         l2g_indices[partial_sum + s]  = n_global_lambda + s;
405674ae819SStefano Zampini         cols_B_delta[partial_sum + s] = dual_dofs_boundary_indices[i];
406674ae819SStefano Zampini         vals_B_delta[partial_sum + s] = 0.0;
407674ae819SStefano Zampini       }
408674ae819SStefano Zampini       /* B_delta */
409674ae819SStefano Zampini       if (n_neg_values > 0) { /* there's a rank next to me to the left */
410674ae819SStefano Zampini         vals_B_delta[partial_sum + n_neg_values - 1] = -1.0;
411674ae819SStefano Zampini       }
412674ae819SStefano Zampini       if (n_neg_values < j) { /* there's a rank next to me to the right */
413674ae819SStefano Zampini         vals_B_delta[partial_sum + n_neg_values] = 1.0;
414674ae819SStefano Zampini       }
415674ae819SStefano Zampini       /* scaling as in Klawonn-Widlund 1999 */
41625d06dbeSstefano_zampini       if (!pcbddc->use_deluxe_scaling) {
417674ae819SStefano Zampini         for (s = 0; s < n_neg_values; s++) {
418674ae819SStefano Zampini           scalar_value = 0.0;
41925d06dbeSstefano_zampini           for (k = 0; k < s + 1; k++) scalar_value += array[k];
420674ae819SStefano Zampini           scaling_factors[partial_sum + s] = -scalar_value;
421674ae819SStefano Zampini         }
422674ae819SStefano Zampini         for (s = 0; s < n_pos_values; s++) {
423674ae819SStefano Zampini           scalar_value = 0.0;
42425d06dbeSstefano_zampini           for (k = s + n_neg_values; k < j; k++) scalar_value += array[k];
425674ae819SStefano Zampini           scaling_factors[partial_sum + s + n_neg_values] = scalar_value;
426674ae819SStefano Zampini         }
42725d06dbeSstefano_zampini       }
428674ae819SStefano Zampini       partial_sum += j;
429674ae819SStefano Zampini     }
430dc456d91SStefano Zampini     cum += aux_local_numbering_2[i];
431674ae819SStefano Zampini   }
4329566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(subset_n, &aux_global_numbering));
4339566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&subset_mult));
4349566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&subset_n));
4359566063dSJacob Faibussowitsch   PetscCall(PetscFree(aux_sums));
4369566063dSJacob Faibussowitsch   PetscCall(PetscFree(aux_local_numbering_1));
4379566063dSJacob Faibussowitsch   PetscCall(PetscFree(dual_dofs_boundary_indices));
43825d06dbeSstefano_zampini   if (all_factors) {
4399566063dSJacob Faibussowitsch     PetscCall(PetscFree(all_factors[0]));
4409566063dSJacob Faibussowitsch     PetscCall(PetscFree(all_factors));
44125d06dbeSstefano_zampini   }
442*9de2952eSStefano Zampini   if (pcis->n) PetscCall(PetscFree(neighbours_set[0]));
443*9de2952eSStefano Zampini   PetscCall(PetscFree2(count, neighbours_set));
444674ae819SStefano Zampini 
445674ae819SStefano Zampini   /* Create local part of B_delta */
4469566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_SELF, &fetidpmat_ctx->B_delta));
4479566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(fetidpmat_ctx->B_delta, n_local_lambda, pcis->n_B, n_local_lambda, pcis->n_B));
4489566063dSJacob Faibussowitsch   PetscCall(MatSetType(fetidpmat_ctx->B_delta, MATSEQAIJ));
4499566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(fetidpmat_ctx->B_delta, 1, NULL));
4509566063dSJacob Faibussowitsch   PetscCall(MatSetOption(fetidpmat_ctx->B_delta, MAT_IGNORE_ZERO_ENTRIES, PETSC_TRUE));
45148a46eb9SPierre Jolivet   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));
4529566063dSJacob Faibussowitsch   PetscCall(PetscFree(vals_B_delta));
4539566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(fetidpmat_ctx->B_delta, MAT_FINAL_ASSEMBLY));
4549566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(fetidpmat_ctx->B_delta, MAT_FINAL_ASSEMBLY));
455674ae819SStefano Zampini 
45625d06dbeSstefano_zampini   BD1 = NULL;
45725d06dbeSstefano_zampini   BD2 = NULL;
458674ae819SStefano Zampini   if (fully_redundant) {
45928b400f6SJacob Faibussowitsch     PetscCheck(!pcbddc->use_deluxe_scaling, comm, PETSC_ERR_SUP, "Deluxe FETIDP with fully-redundant multipliers to be implemented");
4609566063dSJacob Faibussowitsch     PetscCall(MatCreate(PETSC_COMM_SELF, &ScalingMat));
4619566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(ScalingMat, n_local_lambda, n_local_lambda, n_local_lambda, n_local_lambda));
4629566063dSJacob Faibussowitsch     PetscCall(MatSetType(ScalingMat, MATSEQAIJ));
4639566063dSJacob Faibussowitsch     PetscCall(MatSeqAIJSetPreallocation(ScalingMat, 1, NULL));
46448a46eb9SPierre Jolivet     for (i = 0; i < n_local_lambda; i++) PetscCall(MatSetValue(ScalingMat, i, i, scaling_factors[i], INSERT_VALUES));
4659566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(ScalingMat, MAT_FINAL_ASSEMBLY));
4669566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(ScalingMat, MAT_FINAL_ASSEMBLY));
4679566063dSJacob Faibussowitsch     PetscCall(MatMatMult(ScalingMat, fetidpmat_ctx->B_delta, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &fetidpmat_ctx->B_Ddelta));
4689566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&ScalingMat));
469674ae819SStefano Zampini   } else {
4709566063dSJacob Faibussowitsch     PetscCall(MatCreate(PETSC_COMM_SELF, &fetidpmat_ctx->B_Ddelta));
4719566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(fetidpmat_ctx->B_Ddelta, n_local_lambda, pcis->n_B, n_local_lambda, pcis->n_B));
472524221abSStefano Zampini     if (!pcbddc->use_deluxe_scaling || !pcbddc->sub_schurs) {
4739566063dSJacob Faibussowitsch       PetscCall(MatSetType(fetidpmat_ctx->B_Ddelta, MATSEQAIJ));
4749566063dSJacob Faibussowitsch       PetscCall(MatSeqAIJSetPreallocation(fetidpmat_ctx->B_Ddelta, 1, NULL));
47548a46eb9SPierre Jolivet       for (i = 0; i < n_local_lambda; i++) PetscCall(MatSetValue(fetidpmat_ctx->B_Ddelta, i, cols_B_delta[i], scaling_factors[i], INSERT_VALUES));
4769566063dSJacob Faibussowitsch       PetscCall(MatAssemblyBegin(fetidpmat_ctx->B_Ddelta, MAT_FINAL_ASSEMBLY));
4779566063dSJacob Faibussowitsch       PetscCall(MatAssemblyEnd(fetidpmat_ctx->B_Ddelta, MAT_FINAL_ASSEMBLY));
47825d06dbeSstefano_zampini     } else {
47925d06dbeSstefano_zampini       /* scaling as in Klawonn-Widlund 1999 */
48025d06dbeSstefano_zampini       PCBDDCDeluxeScaling deluxe_ctx = pcbddc->deluxe_ctx;
48125d06dbeSstefano_zampini       PCBDDCSubSchurs     sub_schurs = pcbddc->sub_schurs;
48225d06dbeSstefano_zampini       Mat                 T;
48325d06dbeSstefano_zampini       PetscScalar        *W, lwork, *Bwork;
484451a39c7SStefano Zampini       const PetscInt     *idxs = NULL;
48525d06dbeSstefano_zampini       PetscInt            cum, mss, *nnz;
48625d06dbeSstefano_zampini       PetscBLASInt       *pivots, B_lwork, B_N, B_ierr;
48725d06dbeSstefano_zampini 
48828b400f6SJacob Faibussowitsch       PetscCheck(pcbddc->deluxe_singlemat, comm, PETSC_ERR_USER, "Cannot compute B_Ddelta! rerun with -pc_bddc_deluxe_singlemat");
48925d06dbeSstefano_zampini       mss = 0;
4909566063dSJacob Faibussowitsch       PetscCall(PetscCalloc1(pcis->n_B, &nnz));
49125d06dbeSstefano_zampini       if (sub_schurs->is_Ej_all) {
4929566063dSJacob Faibussowitsch         PetscCall(ISGetIndices(sub_schurs->is_Ej_all, &idxs));
49325d06dbeSstefano_zampini         for (i = 0, cum = 0; i < sub_schurs->n_subs; i++) {
49425d06dbeSstefano_zampini           PetscInt subset_size;
49525d06dbeSstefano_zampini 
4969566063dSJacob Faibussowitsch           PetscCall(ISGetLocalSize(sub_schurs->is_subs[i], &subset_size));
49725d06dbeSstefano_zampini           for (j = 0; j < subset_size; j++) nnz[idxs[j + cum]] = subset_size;
49825d06dbeSstefano_zampini           mss = PetscMax(mss, subset_size);
49925d06dbeSstefano_zampini           cum += subset_size;
50025d06dbeSstefano_zampini         }
50125d06dbeSstefano_zampini       }
5029566063dSJacob Faibussowitsch       PetscCall(MatCreate(PETSC_COMM_SELF, &T));
5039566063dSJacob Faibussowitsch       PetscCall(MatSetSizes(T, pcis->n_B, pcis->n_B, pcis->n_B, pcis->n_B));
5049566063dSJacob Faibussowitsch       PetscCall(MatSetType(T, MATSEQAIJ));
5059566063dSJacob Faibussowitsch       PetscCall(MatSeqAIJSetPreallocation(T, 0, nnz));
5069566063dSJacob Faibussowitsch       PetscCall(PetscFree(nnz));
50725d06dbeSstefano_zampini 
50825d06dbeSstefano_zampini       /* workspace allocation */
5097f00194aSstefano_zampini       B_lwork = 0;
5107f00194aSstefano_zampini       if (mss) {
511be8bcbd0Sstefano_zampini         PetscScalar dummy = 1;
512be8bcbd0Sstefano_zampini 
51325d06dbeSstefano_zampini         B_lwork = -1;
5149566063dSJacob Faibussowitsch         PetscCall(PetscBLASIntCast(mss, &B_N));
5159566063dSJacob Faibussowitsch         PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
516792fecdfSBarry Smith         PetscCallBLAS("LAPACKgetri", LAPACKgetri_(&B_N, &dummy, &B_N, &B_N, &lwork, &B_lwork, &B_ierr));
5179566063dSJacob Faibussowitsch         PetscCall(PetscFPTrapPop());
51828b400f6SJacob Faibussowitsch         PetscCheck(!B_ierr, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in query to GETRI Lapack routine %d", (int)B_ierr);
5199566063dSJacob Faibussowitsch         PetscCall(PetscBLASIntCast((PetscInt)PetscRealPart(lwork), &B_lwork));
5207f00194aSstefano_zampini       }
5219566063dSJacob Faibussowitsch       PetscCall(PetscMalloc3(mss * mss, &W, mss, &pivots, B_lwork, &Bwork));
52225d06dbeSstefano_zampini 
52325d06dbeSstefano_zampini       for (i = 0, cum = 0; i < sub_schurs->n_subs; i++) {
5241683a169SBarry Smith         const PetscScalar *M;
52525d06dbeSstefano_zampini         PetscInt           subset_size;
52625d06dbeSstefano_zampini 
5279566063dSJacob Faibussowitsch         PetscCall(ISGetLocalSize(sub_schurs->is_subs[i], &subset_size));
5289566063dSJacob Faibussowitsch         PetscCall(PetscBLASIntCast(subset_size, &B_N));
5299566063dSJacob Faibussowitsch         PetscCall(MatDenseGetArrayRead(deluxe_ctx->seq_mat[i], &M));
5309566063dSJacob Faibussowitsch         PetscCall(PetscArraycpy(W, M, subset_size * subset_size));
5319566063dSJacob Faibussowitsch         PetscCall(MatDenseRestoreArrayRead(deluxe_ctx->seq_mat[i], &M));
5329566063dSJacob Faibussowitsch         PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
533792fecdfSBarry Smith         PetscCallBLAS("LAPACKgetrf", LAPACKgetrf_(&B_N, &B_N, W, &B_N, pivots, &B_ierr));
53428b400f6SJacob Faibussowitsch         PetscCheck(!B_ierr, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in GETRF Lapack routine %d", (int)B_ierr);
535792fecdfSBarry Smith         PetscCallBLAS("LAPACKgetri", LAPACKgetri_(&B_N, W, &B_N, pivots, Bwork, &B_lwork, &B_ierr));
53628b400f6SJacob Faibussowitsch         PetscCheck(!B_ierr, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in GETRI Lapack routine %d", (int)B_ierr);
5379566063dSJacob Faibussowitsch         PetscCall(PetscFPTrapPop());
538451a39c7SStefano Zampini         /* silent static analyzer */
53928b400f6SJacob Faibussowitsch         PetscCheck(idxs, PETSC_COMM_SELF, PETSC_ERR_PLIB, "IDXS not present");
5409566063dSJacob Faibussowitsch         PetscCall(MatSetValues(T, subset_size, idxs + cum, subset_size, idxs + cum, W, INSERT_VALUES));
54125d06dbeSstefano_zampini         cum += subset_size;
54225d06dbeSstefano_zampini       }
5439566063dSJacob Faibussowitsch       PetscCall(MatAssemblyBegin(T, MAT_FINAL_ASSEMBLY));
5449566063dSJacob Faibussowitsch       PetscCall(MatAssemblyEnd(T, MAT_FINAL_ASSEMBLY));
5459566063dSJacob Faibussowitsch       PetscCall(MatMatTransposeMult(T, fetidpmat_ctx->B_delta, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &BD1));
5469566063dSJacob Faibussowitsch       PetscCall(MatMatMult(fetidpmat_ctx->B_delta, BD1, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &BD2));
5479566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&T));
5489566063dSJacob Faibussowitsch       PetscCall(PetscFree3(W, pivots, Bwork));
54948a46eb9SPierre Jolivet       if (sub_schurs->is_Ej_all) PetscCall(ISRestoreIndices(sub_schurs->is_Ej_all, &idxs));
55025d06dbeSstefano_zampini     }
551674ae819SStefano Zampini   }
5529566063dSJacob Faibussowitsch   PetscCall(PetscFree(scaling_factors));
5539566063dSJacob Faibussowitsch   PetscCall(PetscFree(cols_B_delta));
554674ae819SStefano Zampini 
555457d4a33Sstefano_zampini   /* Layout of multipliers */
5569566063dSJacob Faibussowitsch   PetscCall(PetscLayoutCreate(comm, &llay));
5579566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetBlockSize(llay, 1));
5589566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetSize(llay, fetidpmat_ctx->n_lambda));
5599566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(llay));
5609566063dSJacob Faibussowitsch   PetscCall(PetscLayoutGetLocalSize(llay, &fetidpmat_ctx->n));
561457d4a33Sstefano_zampini 
562457d4a33Sstefano_zampini   /* Local work vector of multipliers */
5639566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_SELF, &fetidpmat_ctx->lambda_local));
5649566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(fetidpmat_ctx->lambda_local, n_local_lambda, n_local_lambda));
5659566063dSJacob Faibussowitsch   PetscCall(VecSetType(fetidpmat_ctx->lambda_local, VECSEQ));
566457d4a33Sstefano_zampini 
56725d06dbeSstefano_zampini   if (BD2) {
56825d06dbeSstefano_zampini     ISLocalToGlobalMapping l2g;
56925d06dbeSstefano_zampini     Mat                    T, TA, *pT;
57025d06dbeSstefano_zampini     IS                     is;
57125d06dbeSstefano_zampini     PetscInt               nl, N;
57225d06dbeSstefano_zampini     BDdelta_DN             ctx;
57325d06dbeSstefano_zampini 
5749566063dSJacob Faibussowitsch     PetscCall(PetscLayoutGetLocalSize(llay, &nl));
5759566063dSJacob Faibussowitsch     PetscCall(PetscLayoutGetSize(llay, &N));
5769566063dSJacob Faibussowitsch     PetscCall(MatCreate(comm, &T));
5779566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(T, nl, nl, N, N));
5789566063dSJacob Faibussowitsch     PetscCall(MatSetType(T, MATIS));
5799566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingCreate(comm, 1, n_local_lambda, l2g_indices, PETSC_COPY_VALUES, &l2g));
5809566063dSJacob Faibussowitsch     PetscCall(MatSetLocalToGlobalMapping(T, l2g, l2g));
5819566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingDestroy(&l2g));
5829566063dSJacob Faibussowitsch     PetscCall(MatISSetLocalMat(T, BD2));
5839566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(T, MAT_FINAL_ASSEMBLY));
5849566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(T, MAT_FINAL_ASSEMBLY));
5859566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&BD2));
5869566063dSJacob Faibussowitsch     PetscCall(MatConvert(T, MATAIJ, MAT_INITIAL_MATRIX, &TA));
5879566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&T));
5889566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(comm, n_local_lambda, l2g_indices, PETSC_USE_POINTER, &is));
5899566063dSJacob Faibussowitsch     PetscCall(MatCreateSubMatrices(TA, 1, &is, &is, MAT_INITIAL_MATRIX, &pT));
5909566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&TA));
5919566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&is));
59225d06dbeSstefano_zampini     BD2 = pT[0];
5939566063dSJacob Faibussowitsch     PetscCall(PetscFree(pT));
59425d06dbeSstefano_zampini 
59525d06dbeSstefano_zampini     /* B_Ddelta for non-redundant multipliers with deluxe scaling */
5969566063dSJacob Faibussowitsch     PetscCall(PetscNew(&ctx));
5979566063dSJacob Faibussowitsch     PetscCall(MatSetType(fetidpmat_ctx->B_Ddelta, MATSHELL));
5989566063dSJacob Faibussowitsch     PetscCall(MatShellSetContext(fetidpmat_ctx->B_Ddelta, ctx));
5999566063dSJacob Faibussowitsch     PetscCall(MatShellSetOperation(fetidpmat_ctx->B_Ddelta, MATOP_MULT, (void (*)(void))MatMult_BDdelta_deluxe_nonred));
6009566063dSJacob Faibussowitsch     PetscCall(MatShellSetOperation(fetidpmat_ctx->B_Ddelta, MATOP_MULT_TRANSPOSE, (void (*)(void))MatMultTranspose_BDdelta_deluxe_nonred));
6019566063dSJacob Faibussowitsch     PetscCall(MatShellSetOperation(fetidpmat_ctx->B_Ddelta, MATOP_DESTROY, (void (*)(void))MatDestroy_BDdelta_deluxe_nonred));
6029566063dSJacob Faibussowitsch     PetscCall(MatSetUp(fetidpmat_ctx->B_Ddelta));
60325d06dbeSstefano_zampini 
6049566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)BD1));
60525d06dbeSstefano_zampini     ctx->BD = BD1;
6069566063dSJacob Faibussowitsch     PetscCall(KSPCreate(PETSC_COMM_SELF, &ctx->kBD));
6073821be0aSBarry Smith     PetscCall(KSPSetNestLevel(ctx->kBD, fetidpmat_ctx->pc->kspnestlevel));
6089566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(ctx->kBD, BD2, BD2));
6099566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(fetidpmat_ctx->lambda_local, &ctx->work));
61025d06dbeSstefano_zampini     fetidpmat_ctx->deluxe_nonred = PETSC_TRUE;
61125d06dbeSstefano_zampini   }
6129566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&BD1));
6139566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&BD2));
61425d06dbeSstefano_zampini 
61525d06dbeSstefano_zampini   /* fetidpmat sizes */
61625d06dbeSstefano_zampini   fetidpmat_ctx->n += nPgl;
61725d06dbeSstefano_zampini   fetidpmat_ctx->N = fetidpmat_ctx->n_lambda + nPg;
61825d06dbeSstefano_zampini 
619457d4a33Sstefano_zampini   /* Global vector for FETI-DP linear system */
6209566063dSJacob Faibussowitsch   PetscCall(VecCreate(comm, &fetidp_global));
6219566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(fetidp_global, fetidpmat_ctx->n, fetidpmat_ctx->N));
6229566063dSJacob Faibussowitsch   PetscCall(VecSetType(fetidp_global, VECMPI));
6239566063dSJacob Faibussowitsch   PetscCall(VecSetUp(fetidp_global));
624457d4a33Sstefano_zampini 
6259eec4de8Sstefano_zampini   /* Decide layout for fetidp dofs: if it is a saddle point problem
6269eec4de8Sstefano_zampini      pressure is ordered first in the local part of the global vector
627457d4a33Sstefano_zampini      of the FETI-DP linear system */
628457d4a33Sstefano_zampini   if (nPg) {
629e1214c54Sstefano_zampini     Vec             v;
630af140850Sstefano_zampini     IS              IS_l2g_p, ais;
631457d4a33Sstefano_zampini     PetscLayout     alay;
632457d4a33Sstefano_zampini     const PetscInt *idxs, *pranges, *aranges, *lranges;
633af140850Sstefano_zampini     PetscInt       *l2g_indices_p, rst;
634e1214c54Sstefano_zampini     PetscMPIInt     rank;
635457d4a33Sstefano_zampini 
6369566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nPl, &l2g_indices_p));
6379566063dSJacob Faibussowitsch     PetscCall(VecGetLayout(fetidp_global, &alay));
6389566063dSJacob Faibussowitsch     PetscCall(PetscLayoutGetRanges(alay, &aranges));
6399566063dSJacob Faibussowitsch     PetscCall(PetscLayoutGetRanges(play, &pranges));
6409566063dSJacob Faibussowitsch     PetscCall(PetscLayoutGetRanges(llay, &lranges));
641e1214c54Sstefano_zampini 
6429566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)fetidp_global), &rank));
6439566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PetscObjectComm((PetscObject)fetidp_global), pranges[rank + 1] - pranges[rank], aranges[rank], 1, &fetidpmat_ctx->pressure));
6449566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)fetidpmat_ctx->pressure, "F_P"));
6459566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PetscObjectComm((PetscObject)fetidp_global), lranges[rank + 1] - lranges[rank], aranges[rank] + pranges[rank + 1] - pranges[rank], 1, &fetidpmat_ctx->lagrange));
6469566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)fetidpmat_ctx->lagrange, "F_L"));
6479566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingGetIndices(l2gmap_p, &idxs));
648af140850Sstefano_zampini     /* shift local to global indices for pressure */
649457d4a33Sstefano_zampini     for (i = 0; i < nPl; i++) {
650131c27b5Sprj-       PetscMPIInt owner;
651457d4a33Sstefano_zampini 
6529566063dSJacob Faibussowitsch       PetscCall(PetscLayoutFindOwner(play, idxs[i], &owner));
653457d4a33Sstefano_zampini       l2g_indices_p[i] = idxs[i] - pranges[owner] + aranges[owner];
654457d4a33Sstefano_zampini     }
6559566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingRestoreIndices(l2gmap_p, &idxs));
6569566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(comm, nPl, l2g_indices_p, PETSC_OWN_POINTER, &IS_l2g_p));
657af140850Sstefano_zampini 
658e1214c54Sstefano_zampini     /* local to global scatter for pressure */
6599566063dSJacob Faibussowitsch     PetscCall(VecScatterCreate(fetidpmat_ctx->vP, NULL, fetidp_global, IS_l2g_p, &fetidpmat_ctx->l2g_p));
6609566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&IS_l2g_p));
661457d4a33Sstefano_zampini 
662e1214c54Sstefano_zampini     /* scatter for lagrange multipliers only */
6639566063dSJacob Faibussowitsch     PetscCall(VecCreate(comm, &v));
6649566063dSJacob Faibussowitsch     PetscCall(VecSetType(v, VECSTANDARD));
6659566063dSJacob Faibussowitsch     PetscCall(VecSetLayout(v, llay));
6669566063dSJacob Faibussowitsch     PetscCall(VecSetUp(v));
6679566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(comm, n_local_lambda, l2g_indices, PETSC_COPY_VALUES, &ais));
6689566063dSJacob Faibussowitsch     PetscCall(VecScatterCreate(fetidpmat_ctx->lambda_local, NULL, v, ais, &fetidpmat_ctx->l2g_lambda_only));
6699566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&ais));
6709566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&v));
671e1214c54Sstefano_zampini 
672af140850Sstefano_zampini     /* shift local to global indices for multipliers */
673457d4a33Sstefano_zampini     for (i = 0; i < n_local_lambda; i++) {
674131c27b5Sprj-       PetscInt    ps;
675131c27b5Sprj-       PetscMPIInt owner;
676457d4a33Sstefano_zampini 
6779566063dSJacob Faibussowitsch       PetscCall(PetscLayoutFindOwner(llay, l2g_indices[i], &owner));
678457d4a33Sstefano_zampini       ps             = pranges[owner + 1] - pranges[owner];
679457d4a33Sstefano_zampini       l2g_indices[i] = l2g_indices[i] - lranges[owner] + aranges[owner] + ps;
680457d4a33Sstefano_zampini     }
681457d4a33Sstefano_zampini 
682e1214c54Sstefano_zampini     /* scatter from alldofs to pressures global fetidp vector */
6839566063dSJacob Faibussowitsch     PetscCall(PetscLayoutGetRange(alay, &rst, NULL));
6849566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(comm, nPgl, rst, 1, &ais));
6859566063dSJacob Faibussowitsch     PetscCall(VecScatterCreate(pcis->vec1_global, pP, fetidp_global, ais, &fetidpmat_ctx->g2g_p));
6869566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&ais));
687457d4a33Sstefano_zampini   }
6889566063dSJacob Faibussowitsch   PetscCall(PetscLayoutDestroy(&llay));
6899566063dSJacob Faibussowitsch   PetscCall(PetscLayoutDestroy(&play));
6909566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(comm, n_local_lambda, l2g_indices, PETSC_OWN_POINTER, &IS_l2g_lambda));
691a1c0d0daSstefano_zampini 
6929cc7774eSstefano_zampini   /* scatter from local to global multipliers */
6939566063dSJacob Faibussowitsch   PetscCall(VecScatterCreate(fetidpmat_ctx->lambda_local, NULL, fetidp_global, IS_l2g_lambda, &fetidpmat_ctx->l2g_lambda));
6949566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&IS_l2g_lambda));
6959566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingDestroy(&l2gmap_p));
6969566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&fetidp_global));
697457d4a33Sstefano_zampini 
698a1c0d0daSstefano_zampini   /* Create some work vectors needed by fetidp */
6999566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(pcis->vec1_B, &fetidpmat_ctx->temp_solution_B));
7009566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(pcis->vec1_D, &fetidpmat_ctx->temp_solution_D));
7013ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
702674ae819SStefano Zampini }
703674ae819SStefano Zampini 
704d71ae5a4SJacob Faibussowitsch PetscErrorCode PCBDDCSetupFETIDPPCContext(Mat fetimat, FETIDPPC_ctx fetidppc_ctx)
705d71ae5a4SJacob Faibussowitsch {
706674ae819SStefano Zampini   FETIDPMat_ctx mat_ctx;
70725d06dbeSstefano_zampini   PC_BDDC      *pcbddc = (PC_BDDC *)fetidppc_ctx->pc->data;
70825d06dbeSstefano_zampini   PC_IS        *pcis   = (PC_IS *)fetidppc_ctx->pc->data;
709f28b6018SStefano Zampini   PetscBool     lumped = PETSC_FALSE;
710674ae819SStefano Zampini 
711674ae819SStefano Zampini   PetscFunctionBegin;
7129566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(fetimat, &mat_ctx));
713674ae819SStefano Zampini   /* get references from objects created when setting up feti mat context */
7149566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)mat_ctx->lambda_local));
715674ae819SStefano Zampini   fetidppc_ctx->lambda_local = mat_ctx->lambda_local;
7169566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)mat_ctx->B_Ddelta));
717674ae819SStefano Zampini   fetidppc_ctx->B_Ddelta = mat_ctx->B_Ddelta;
71825d06dbeSstefano_zampini   if (mat_ctx->deluxe_nonred) {
71925d06dbeSstefano_zampini     PC            pc, mpc;
72025d06dbeSstefano_zampini     BDdelta_DN    ctx;
7213ca39a21SBarry Smith     MatSolverType solver;
72225d06dbeSstefano_zampini     const char   *prefix;
72325d06dbeSstefano_zampini 
7249566063dSJacob Faibussowitsch     PetscCall(MatShellGetContext(mat_ctx->B_Ddelta, &ctx));
7259566063dSJacob Faibussowitsch     PetscCall(KSPSetType(ctx->kBD, KSPPREONLY));
7269566063dSJacob Faibussowitsch     PetscCall(KSPGetPC(ctx->kBD, &mpc));
7279566063dSJacob Faibussowitsch     PetscCall(KSPGetPC(pcbddc->ksp_D, &pc));
7289566063dSJacob Faibussowitsch     PetscCall(PCSetType(mpc, PCLU));
7299566063dSJacob Faibussowitsch     PetscCall(PCFactorGetMatSolverType(pc, (MatSolverType *)&solver));
7301baa6e33SBarry Smith     if (solver) PetscCall(PCFactorSetMatSolverType(mpc, solver));
7319566063dSJacob Faibussowitsch     PetscCall(MatGetOptionsPrefix(fetimat, &prefix));
7329566063dSJacob Faibussowitsch     PetscCall(KSPSetOptionsPrefix(ctx->kBD, prefix));
7339566063dSJacob Faibussowitsch     PetscCall(KSPAppendOptionsPrefix(ctx->kBD, "bddelta_"));
7349566063dSJacob Faibussowitsch     PetscCall(KSPSetFromOptions(ctx->kBD));
73525d06dbeSstefano_zampini   }
73625d06dbeSstefano_zampini 
737e1214c54Sstefano_zampini   if (mat_ctx->l2g_lambda_only) {
7389566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)mat_ctx->l2g_lambda_only));
739e1214c54Sstefano_zampini     fetidppc_ctx->l2g_lambda = mat_ctx->l2g_lambda_only;
740e1214c54Sstefano_zampini   } else {
7419566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)mat_ctx->l2g_lambda));
742674ae819SStefano Zampini     fetidppc_ctx->l2g_lambda = mat_ctx->l2g_lambda;
743e1214c54Sstefano_zampini   }
744f28b6018SStefano Zampini   /* Dirichlet preconditioner */
7459566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, ((PetscObject)fetimat)->prefix, "-pc_lumped", &lumped, NULL));
746f28b6018SStefano Zampini   if (!lumped) {
7479feb908dSstefano_zampini     IS        iV;
7489c2d02cdSstefano_zampini     PetscBool discrete_harmonic = PETSC_FALSE;
7499c2d02cdSstefano_zampini 
7509566063dSJacob Faibussowitsch     PetscCall(PetscObjectQuery((PetscObject)fetidppc_ctx->pc, "__KSPFETIDP_iV", (PetscObject *)&iV));
75148a46eb9SPierre Jolivet     if (iV) PetscCall(PetscOptionsGetBool(NULL, ((PetscObject)fetimat)->prefix, "-pc_discrete_harmonic", &discrete_harmonic, NULL));
7529c2d02cdSstefano_zampini     if (discrete_harmonic) {
7539c2d02cdSstefano_zampini       KSP             sksp;
7549c2d02cdSstefano_zampini       PC              pc;
755111315fdSstefano_zampini       PCBDDCSubSchurs sub_schurs = pcbddc->sub_schurs;
7569c2d02cdSstefano_zampini       Mat             A_II, A_IB, A_BI;
757111315fdSstefano_zampini       IS              iP = NULL;
758111315fdSstefano_zampini       PetscBool       isshell, reuse = PETSC_FALSE;
7599c2d02cdSstefano_zampini       KSPType         ksptype;
7609c2d02cdSstefano_zampini       const char     *prefix;
7619c2d02cdSstefano_zampini 
7629c2d02cdSstefano_zampini       /*
7639c2d02cdSstefano_zampini         We constructs a Schur complement for
7649c2d02cdSstefano_zampini 
7659c2d02cdSstefano_zampini         | A_II A_ID |
7669c2d02cdSstefano_zampini         | A_DI A_DD |
7679c2d02cdSstefano_zampini 
7689c2d02cdSstefano_zampini         instead of
7699c2d02cdSstefano_zampini 
7709c2d02cdSstefano_zampini         | A_II  B^t_II A_ID |
7719c2d02cdSstefano_zampini         | B_II -C_II   B_ID |
7729c2d02cdSstefano_zampini         | A_DI  B^t_ID A_DD |
7739c2d02cdSstefano_zampini 
7749c2d02cdSstefano_zampini       */
775111315fdSstefano_zampini       if (sub_schurs && sub_schurs->reuse_solver) {
7769566063dSJacob Faibussowitsch         PetscCall(PetscObjectQuery((PetscObject)sub_schurs->A, "__KSPFETIDP_iP", (PetscObject *)&iP));
777111315fdSstefano_zampini         if (iP) reuse = PETSC_TRUE;
778111315fdSstefano_zampini       }
779111315fdSstefano_zampini       if (!reuse) {
780111315fdSstefano_zampini         IS       aB;
781111315fdSstefano_zampini         PetscInt nb;
7829566063dSJacob Faibussowitsch         PetscCall(ISGetLocalSize(pcis->is_B_local, &nb));
7839566063dSJacob Faibussowitsch         PetscCall(ISCreateStride(PetscObjectComm((PetscObject)pcis->A_II), nb, 0, 1, &aB));
7849566063dSJacob Faibussowitsch         PetscCall(MatCreateSubMatrix(pcis->A_II, iV, iV, MAT_INITIAL_MATRIX, &A_II));
7859566063dSJacob Faibussowitsch         PetscCall(MatCreateSubMatrix(pcis->A_IB, iV, aB, MAT_INITIAL_MATRIX, &A_IB));
7869566063dSJacob Faibussowitsch         PetscCall(MatCreateSubMatrix(pcis->A_BI, aB, iV, MAT_INITIAL_MATRIX, &A_BI));
7879566063dSJacob Faibussowitsch         PetscCall(ISDestroy(&aB));
788111315fdSstefano_zampini       } else {
7899566063dSJacob Faibussowitsch         PetscCall(MatCreateSubMatrix(sub_schurs->A, pcis->is_I_local, pcis->is_B_local, MAT_INITIAL_MATRIX, &A_IB));
7909566063dSJacob Faibussowitsch         PetscCall(MatCreateSubMatrix(sub_schurs->A, pcis->is_B_local, pcis->is_I_local, MAT_INITIAL_MATRIX, &A_BI));
7919566063dSJacob Faibussowitsch         PetscCall(PetscObjectReference((PetscObject)pcis->A_II));
792111315fdSstefano_zampini         A_II = pcis->A_II;
793111315fdSstefano_zampini       }
7949566063dSJacob Faibussowitsch       PetscCall(MatCreateSchurComplement(A_II, A_II, A_IB, A_BI, pcis->A_BB, &fetidppc_ctx->S_j));
7959c2d02cdSstefano_zampini 
7969c2d02cdSstefano_zampini       /* propagate settings of solver */
7979566063dSJacob Faibussowitsch       PetscCall(MatSchurComplementGetKSP(fetidppc_ctx->S_j, &sksp));
7989566063dSJacob Faibussowitsch       PetscCall(KSPGetType(pcis->ksp_D, &ksptype));
7999566063dSJacob Faibussowitsch       PetscCall(KSPSetType(sksp, ksptype));
8009566063dSJacob Faibussowitsch       PetscCall(KSPGetPC(pcis->ksp_D, &pc));
8019566063dSJacob Faibussowitsch       PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCSHELL, &isshell));
8025334bea6Sstefano_zampini       if (!isshell) {
8033ca39a21SBarry Smith         MatSolverType solver;
8045334bea6Sstefano_zampini         PCType        pctype;
8055334bea6Sstefano_zampini 
8069566063dSJacob Faibussowitsch         PetscCall(PCGetType(pc, &pctype));
8079566063dSJacob Faibussowitsch         PetscCall(PCFactorGetMatSolverType(pc, (MatSolverType *)&solver));
8089566063dSJacob Faibussowitsch         PetscCall(KSPGetPC(sksp, &pc));
8099566063dSJacob Faibussowitsch         PetscCall(PCSetType(pc, pctype));
8101baa6e33SBarry Smith         if (solver) PetscCall(PCFactorSetMatSolverType(pc, solver));
8115334bea6Sstefano_zampini       } else {
8129566063dSJacob Faibussowitsch         PetscCall(KSPGetPC(sksp, &pc));
8139566063dSJacob Faibussowitsch         PetscCall(PCSetType(pc, PCLU));
8145334bea6Sstefano_zampini       }
8159566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&A_II));
8169566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&A_IB));
8179566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&A_BI));
8189566063dSJacob Faibussowitsch       PetscCall(MatGetOptionsPrefix(fetimat, &prefix));
8199566063dSJacob Faibussowitsch       PetscCall(KSPSetOptionsPrefix(sksp, prefix));
8209566063dSJacob Faibussowitsch       PetscCall(KSPAppendOptionsPrefix(sksp, "harmonic_"));
8219566063dSJacob Faibussowitsch       PetscCall(KSPSetFromOptions(sksp));
822111315fdSstefano_zampini       if (reuse) {
8239566063dSJacob Faibussowitsch         PetscCall(KSPSetPC(sksp, sub_schurs->reuse_solver->interior_solver));
8249566063dSJacob Faibussowitsch         PetscCall(PetscObjectIncrementTabLevel((PetscObject)sub_schurs->reuse_solver->interior_solver, (PetscObject)sksp, 0));
825111315fdSstefano_zampini       }
8269c2d02cdSstefano_zampini     } else { /* default Dirichlet preconditioner is pde-harmonic */
8279566063dSJacob Faibussowitsch       PetscCall(MatCreateSchurComplement(pcis->A_II, pcis->A_II, pcis->A_IB, pcis->A_BI, pcis->A_BB, &fetidppc_ctx->S_j));
8289566063dSJacob Faibussowitsch       PetscCall(MatSchurComplementSetKSP(fetidppc_ctx->S_j, pcis->ksp_D));
8299c2d02cdSstefano_zampini     }
830f28b6018SStefano Zampini   } else {
8319566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)pcis->A_BB));
832f28b6018SStefano Zampini     fetidppc_ctx->S_j = pcis->A_BB;
833f28b6018SStefano Zampini   }
834af140850Sstefano_zampini   /* saddle-point */
835af140850Sstefano_zampini   if (mat_ctx->xPg) {
8369566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)mat_ctx->xPg));
837af140850Sstefano_zampini     fetidppc_ctx->xPg = mat_ctx->xPg;
8389566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)mat_ctx->yPg));
839af140850Sstefano_zampini     fetidppc_ctx->yPg = mat_ctx->yPg;
840af140850Sstefano_zampini   }
8413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
842674ae819SStefano Zampini }
843674ae819SStefano Zampini 
84466976f2fSJacob Faibussowitsch static PetscErrorCode FETIDPMatMult_Kernel(Mat fetimat, Vec x, Vec y, PetscBool trans)
845d71ae5a4SJacob Faibussowitsch {
846674ae819SStefano Zampini   FETIDPMat_ctx mat_ctx;
847617d11aeSStefano Zampini   PC_BDDC      *pcbddc;
848674ae819SStefano Zampini   PC_IS        *pcis;
849674ae819SStefano Zampini 
850674ae819SStefano Zampini   PetscFunctionBegin;
8519566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(fetimat, &mat_ctx));
852674ae819SStefano Zampini   pcis   = (PC_IS *)mat_ctx->pc->data;
853617d11aeSStefano Zampini   pcbddc = (PC_BDDC *)mat_ctx->pc->data;
854674ae819SStefano Zampini   /* Application of B_delta^T */
8559566063dSJacob Faibussowitsch   PetscCall(VecSet(pcis->vec1_B, 0.));
8569566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(mat_ctx->l2g_lambda, x, mat_ctx->lambda_local, INSERT_VALUES, SCATTER_REVERSE));
8579566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(mat_ctx->l2g_lambda, x, mat_ctx->lambda_local, INSERT_VALUES, SCATTER_REVERSE));
8589566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(mat_ctx->B_delta, mat_ctx->lambda_local, pcis->vec1_B));
859af140850Sstefano_zampini 
860af140850Sstefano_zampini   /* Add contribution from saddle point */
861af140850Sstefano_zampini   if (mat_ctx->l2g_p) {
8629566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(mat_ctx->l2g_p, x, mat_ctx->vP, INSERT_VALUES, SCATTER_REVERSE));
8639566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(mat_ctx->l2g_p, x, mat_ctx->vP, INSERT_VALUES, SCATTER_REVERSE));
864af140850Sstefano_zampini     if (pcbddc->switch_static) {
865ef425aeaSstefano_zampini       if (trans) {
8669566063dSJacob Faibussowitsch         PetscCall(MatMultTranspose(mat_ctx->B_BI, mat_ctx->vP, pcis->vec1_D));
867ef425aeaSstefano_zampini       } else {
8689566063dSJacob Faibussowitsch         PetscCall(MatMult(mat_ctx->Bt_BI, mat_ctx->vP, pcis->vec1_D));
869af140850Sstefano_zampini       }
870ef425aeaSstefano_zampini     }
871ef425aeaSstefano_zampini     if (trans) {
8729566063dSJacob Faibussowitsch       PetscCall(MatMultTransposeAdd(mat_ctx->B_BB, mat_ctx->vP, pcis->vec1_B, pcis->vec1_B));
873ef425aeaSstefano_zampini     } else {
8749566063dSJacob Faibussowitsch       PetscCall(MatMultAdd(mat_ctx->Bt_BB, mat_ctx->vP, pcis->vec1_B, pcis->vec1_B));
875ef425aeaSstefano_zampini     }
876af140850Sstefano_zampini   } else {
8771baa6e33SBarry Smith     if (pcbddc->switch_static) PetscCall(VecSet(pcis->vec1_D, 0.0));
878af140850Sstefano_zampini   }
879af140850Sstefano_zampini   /* Application of \widetilde{S}^-1 */
8809566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(pcbddc->benign_p0, pcbddc->benign_n));
8819566063dSJacob Faibussowitsch   PetscCall(PCBDDCApplyInterfacePreconditioner(mat_ctx->pc, trans));
8829566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(pcbddc->benign_p0, pcbddc->benign_n));
8839566063dSJacob Faibussowitsch   PetscCall(VecSet(y, 0.0));
884674ae819SStefano Zampini   /* Application of B_delta */
8859566063dSJacob Faibussowitsch   PetscCall(MatMult(mat_ctx->B_delta, pcis->vec1_B, mat_ctx->lambda_local));
886af140850Sstefano_zampini   /* Contribution from boundary pressures */
887af140850Sstefano_zampini   if (mat_ctx->C) {
888af140850Sstefano_zampini     const PetscScalar *lx;
889af140850Sstefano_zampini     PetscScalar       *ly;
890af140850Sstefano_zampini 
89115579a77SStefano Zampini     /* pressure ordered first in the local part of x and y */
8929566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(x, &lx));
8939566063dSJacob Faibussowitsch     PetscCall(VecGetArray(y, &ly));
8949566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(mat_ctx->xPg, lx));
8959566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(mat_ctx->yPg, ly));
896ef425aeaSstefano_zampini     if (trans) {
8979566063dSJacob Faibussowitsch       PetscCall(MatMultTranspose(mat_ctx->C, mat_ctx->xPg, mat_ctx->yPg));
898ef425aeaSstefano_zampini     } else {
8999566063dSJacob Faibussowitsch       PetscCall(MatMult(mat_ctx->C, mat_ctx->xPg, mat_ctx->yPg));
900ef425aeaSstefano_zampini     }
9019566063dSJacob Faibussowitsch     PetscCall(VecResetArray(mat_ctx->xPg));
9029566063dSJacob Faibussowitsch     PetscCall(VecResetArray(mat_ctx->yPg));
9039566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(x, &lx));
9049566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(y, &ly));
905af140850Sstefano_zampini   }
906af140850Sstefano_zampini   /* Add contribution from saddle point */
907af140850Sstefano_zampini   if (mat_ctx->l2g_p) {
908ef425aeaSstefano_zampini     if (trans) {
9099566063dSJacob Faibussowitsch       PetscCall(MatMultTranspose(mat_ctx->Bt_BB, pcis->vec1_B, mat_ctx->vP));
910ef425aeaSstefano_zampini     } else {
9119566063dSJacob Faibussowitsch       PetscCall(MatMult(mat_ctx->B_BB, pcis->vec1_B, mat_ctx->vP));
912ef425aeaSstefano_zampini     }
913af140850Sstefano_zampini     if (pcbddc->switch_static) {
914ef425aeaSstefano_zampini       if (trans) {
9159566063dSJacob Faibussowitsch         PetscCall(MatMultTransposeAdd(mat_ctx->Bt_BI, pcis->vec1_D, mat_ctx->vP, mat_ctx->vP));
916ef425aeaSstefano_zampini       } else {
9179566063dSJacob Faibussowitsch         PetscCall(MatMultAdd(mat_ctx->B_BI, pcis->vec1_D, mat_ctx->vP, mat_ctx->vP));
918af140850Sstefano_zampini       }
919ef425aeaSstefano_zampini     }
9209566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(mat_ctx->l2g_p, mat_ctx->vP, y, ADD_VALUES, SCATTER_FORWARD));
9219566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(mat_ctx->l2g_p, mat_ctx->vP, y, ADD_VALUES, SCATTER_FORWARD));
922af140850Sstefano_zampini   }
9239566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(mat_ctx->l2g_lambda, mat_ctx->lambda_local, y, ADD_VALUES, SCATTER_FORWARD));
9249566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(mat_ctx->l2g_lambda, mat_ctx->lambda_local, y, ADD_VALUES, SCATTER_FORWARD));
9253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
926674ae819SStefano Zampini }
927674ae819SStefano Zampini 
928d71ae5a4SJacob Faibussowitsch PetscErrorCode FETIDPMatMult(Mat fetimat, Vec x, Vec y)
929d71ae5a4SJacob Faibussowitsch {
930edf7251bSStefano Zampini   PetscFunctionBegin;
9319566063dSJacob Faibussowitsch   PetscCall(FETIDPMatMult_Kernel(fetimat, x, y, PETSC_FALSE));
9323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
933ef425aeaSstefano_zampini }
934ef425aeaSstefano_zampini 
935d71ae5a4SJacob Faibussowitsch PetscErrorCode FETIDPMatMultTranspose(Mat fetimat, Vec x, Vec y)
936d71ae5a4SJacob Faibussowitsch {
937ef425aeaSstefano_zampini   PetscFunctionBegin;
9389566063dSJacob Faibussowitsch   PetscCall(FETIDPMatMult_Kernel(fetimat, x, y, PETSC_TRUE));
9393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
940edf7251bSStefano Zampini }
941edf7251bSStefano Zampini 
94266976f2fSJacob Faibussowitsch static PetscErrorCode FETIDPPCApply_Kernel(PC fetipc, Vec x, Vec y, PetscBool trans)
943d71ae5a4SJacob Faibussowitsch {
944674ae819SStefano Zampini   FETIDPPC_ctx pc_ctx;
945674ae819SStefano Zampini   PC_IS       *pcis;
946674ae819SStefano Zampini 
947674ae819SStefano Zampini   PetscFunctionBegin;
9489566063dSJacob Faibussowitsch   PetscCall(PCShellGetContext(fetipc, &pc_ctx));
949674ae819SStefano Zampini   pcis = (PC_IS *)pc_ctx->pc->data;
950674ae819SStefano Zampini   /* Application of B_Ddelta^T */
9519566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(pc_ctx->l2g_lambda, x, pc_ctx->lambda_local, INSERT_VALUES, SCATTER_REVERSE));
9529566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(pc_ctx->l2g_lambda, x, pc_ctx->lambda_local, INSERT_VALUES, SCATTER_REVERSE));
9539566063dSJacob Faibussowitsch   PetscCall(VecSet(pcis->vec2_B, 0.0));
9549566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(pc_ctx->B_Ddelta, pc_ctx->lambda_local, pcis->vec2_B));
955ed6c3d69SStefano Zampini   /* Application of local Schur complement */
95616597391Sstefano_zampini   if (trans) {
9579566063dSJacob Faibussowitsch     PetscCall(MatMultTranspose(pc_ctx->S_j, pcis->vec2_B, pcis->vec1_B));
95816597391Sstefano_zampini   } else {
9599566063dSJacob Faibussowitsch     PetscCall(MatMult(pc_ctx->S_j, pcis->vec2_B, pcis->vec1_B));
96016597391Sstefano_zampini   }
961edf7251bSStefano Zampini   /* Application of B_Ddelta */
9629566063dSJacob Faibussowitsch   PetscCall(MatMult(pc_ctx->B_Ddelta, pcis->vec1_B, pc_ctx->lambda_local));
9639566063dSJacob Faibussowitsch   PetscCall(VecSet(y, 0.0));
9649566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(pc_ctx->l2g_lambda, pc_ctx->lambda_local, y, ADD_VALUES, SCATTER_FORWARD));
9659566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(pc_ctx->l2g_lambda, pc_ctx->lambda_local, y, ADD_VALUES, SCATTER_FORWARD));
9663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
967edf7251bSStefano Zampini }
968edf7251bSStefano Zampini 
969d71ae5a4SJacob Faibussowitsch PetscErrorCode FETIDPPCApply(PC pc, Vec x, Vec y)
970d71ae5a4SJacob Faibussowitsch {
971edf7251bSStefano Zampini   PetscFunctionBegin;
9729566063dSJacob Faibussowitsch   PetscCall(FETIDPPCApply_Kernel(pc, x, y, PETSC_FALSE));
9733ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
97416597391Sstefano_zampini }
97516597391Sstefano_zampini 
976d71ae5a4SJacob Faibussowitsch PetscErrorCode FETIDPPCApplyTranspose(PC pc, Vec x, Vec y)
977d71ae5a4SJacob Faibussowitsch {
97816597391Sstefano_zampini   PetscFunctionBegin;
9799566063dSJacob Faibussowitsch   PetscCall(FETIDPPCApply_Kernel(pc, x, y, PETSC_TRUE));
9803ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
981674ae819SStefano Zampini }
982c45b8d2dSstefano_zampini 
983d71ae5a4SJacob Faibussowitsch PetscErrorCode FETIDPPCView(PC pc, PetscViewer viewer)
984d71ae5a4SJacob Faibussowitsch {
985c45b8d2dSstefano_zampini   FETIDPPC_ctx pc_ctx;
986c45b8d2dSstefano_zampini   PetscBool    iascii;
987c45b8d2dSstefano_zampini   PetscViewer  sviewer;
988c45b8d2dSstefano_zampini 
989c45b8d2dSstefano_zampini   PetscFunctionBegin;
9909566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
991c45b8d2dSstefano_zampini   if (iascii) {
992c45b8d2dSstefano_zampini     PetscMPIInt rank;
99325d06dbeSstefano_zampini     PetscBool   isschur, isshell;
994c45b8d2dSstefano_zampini 
9959566063dSJacob Faibussowitsch     PetscCall(PCShellGetContext(pc, &pc_ctx));
9969566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank));
9979566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)pc_ctx->S_j, MATSCHURCOMPLEMENT, &isschur));
998c45b8d2dSstefano_zampini     if (isschur) {
9999566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "  Dirichlet preconditioner (just from rank 0)\n"));
1000c45b8d2dSstefano_zampini     } else {
10019566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "  Lumped preconditioner (just from rank 0)\n"));
1002c45b8d2dSstefano_zampini     }
10039566063dSJacob Faibussowitsch     PetscCall(PetscViewerGetSubViewer(viewer, PetscObjectComm((PetscObject)pc_ctx->S_j), &sviewer));
1004dd400576SPatrick Sanan     if (rank == 0) {
10059566063dSJacob Faibussowitsch       PetscCall(PetscViewerPushFormat(sviewer, PETSC_VIEWER_ASCII_INFO));
10069566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPushTab(sviewer));
10079566063dSJacob Faibussowitsch       PetscCall(MatView(pc_ctx->S_j, sviewer));
10089566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPopTab(sviewer));
10099566063dSJacob Faibussowitsch       PetscCall(PetscViewerPopFormat(sviewer));
1010c45b8d2dSstefano_zampini     }
10119566063dSJacob Faibussowitsch     PetscCall(PetscViewerRestoreSubViewer(viewer, PetscObjectComm((PetscObject)pc_ctx->S_j), &sviewer));
10129566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)pc_ctx->B_Ddelta, MATSHELL, &isshell));
101325d06dbeSstefano_zampini     if (isshell) {
101425d06dbeSstefano_zampini       BDdelta_DN ctx;
10159566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "  FETI-DP BDdelta: DB^t * (B D^-1 B^t)^-1 for deluxe scaling (just from rank 0)\n"));
10169566063dSJacob Faibussowitsch       PetscCall(MatShellGetContext(pc_ctx->B_Ddelta, &ctx));
10179566063dSJacob Faibussowitsch       PetscCall(PetscViewerGetSubViewer(viewer, PetscObjectComm((PetscObject)pc_ctx->S_j), &sviewer));
1018dd400576SPatrick Sanan       if (rank == 0) {
101915579a77SStefano Zampini         PetscInt tl;
102015579a77SStefano Zampini 
10219566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIGetTab(sviewer, &tl));
10229566063dSJacob Faibussowitsch         PetscCall(PetscObjectSetTabLevel((PetscObject)ctx->kBD, tl));
10239566063dSJacob Faibussowitsch         PetscCall(KSPView(ctx->kBD, sviewer));
10249566063dSJacob Faibussowitsch         PetscCall(PetscViewerPushFormat(sviewer, PETSC_VIEWER_ASCII_INFO));
10259566063dSJacob Faibussowitsch         PetscCall(MatView(ctx->BD, sviewer));
10269566063dSJacob Faibussowitsch         PetscCall(PetscViewerPopFormat(sviewer));
102725d06dbeSstefano_zampini       }
10289566063dSJacob Faibussowitsch       PetscCall(PetscViewerRestoreSubViewer(viewer, PetscObjectComm((PetscObject)pc_ctx->S_j), &sviewer));
1029e5afcf28SBarry Smith     }
10309566063dSJacob Faibussowitsch     PetscCall(PetscViewerFlush(viewer));
1031c45b8d2dSstefano_zampini   }
10323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1033c45b8d2dSstefano_zampini }
1034