xref: /petsc/src/ksp/pc/impls/bddc/bddcfetidp.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
15e5bbd0aSStefano Zampini #include <petsc/private/pcbddcimpl.h>
25e5bbd0aSStefano Zampini #include <petsc/private/pcbddcprivateimpl.h>
325d06dbeSstefano_zampini #include <petscblaslapack.h>
425d06dbeSstefano_zampini 
59371c9d4SSatish Balay static PetscErrorCode MatMult_BDdelta_deluxe_nonred(Mat A, Vec x, Vec y) {
625d06dbeSstefano_zampini   BDdelta_DN ctx;
725d06dbeSstefano_zampini 
825d06dbeSstefano_zampini   PetscFunctionBegin;
99566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(A, &ctx));
109566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(ctx->BD, x, ctx->work));
119566063dSJacob Faibussowitsch   PetscCall(KSPSolveTranspose(ctx->kBD, ctx->work, y));
12c0decd05SBarry Smith   /* No PC so cannot propagate up failure in KSPSolveTranspose() */
1325d06dbeSstefano_zampini   PetscFunctionReturn(0);
1425d06dbeSstefano_zampini }
1525d06dbeSstefano_zampini 
169371c9d4SSatish Balay static PetscErrorCode MatMultTranspose_BDdelta_deluxe_nonred(Mat A, Vec x, Vec y) {
1725d06dbeSstefano_zampini   BDdelta_DN ctx;
1825d06dbeSstefano_zampini 
1925d06dbeSstefano_zampini   PetscFunctionBegin;
209566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(A, &ctx));
219566063dSJacob Faibussowitsch   PetscCall(KSPSolve(ctx->kBD, x, ctx->work));
22c0decd05SBarry Smith   /* No PC so cannot propagate up failure in KSPSolve() */
239566063dSJacob Faibussowitsch   PetscCall(MatMult(ctx->BD, ctx->work, y));
2425d06dbeSstefano_zampini   PetscFunctionReturn(0);
2525d06dbeSstefano_zampini }
2625d06dbeSstefano_zampini 
279371c9d4SSatish Balay static PetscErrorCode MatDestroy_BDdelta_deluxe_nonred(Mat A) {
2825d06dbeSstefano_zampini   BDdelta_DN ctx;
2925d06dbeSstefano_zampini 
3025d06dbeSstefano_zampini   PetscFunctionBegin;
319566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(A, &ctx));
329566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&ctx->BD));
339566063dSJacob Faibussowitsch   PetscCall(KSPDestroy(&ctx->kBD));
349566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ctx->work));
359566063dSJacob Faibussowitsch   PetscCall(PetscFree(ctx));
3625d06dbeSstefano_zampini   PetscFunctionReturn(0);
3725d06dbeSstefano_zampini }
3825d06dbeSstefano_zampini 
399371c9d4SSatish Balay PetscErrorCode PCBDDCCreateFETIDPMatContext(PC pc, FETIDPMat_ctx *fetidpmat_ctx) {
40674ae819SStefano Zampini   FETIDPMat_ctx newctx;
41674ae819SStefano Zampini 
42674ae819SStefano Zampini   PetscFunctionBegin;
439566063dSJacob Faibussowitsch   PetscCall(PetscNew(&newctx));
44674ae819SStefano Zampini   /* increase the reference count for BDDC preconditioner */
459566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)pc));
46674ae819SStefano Zampini   newctx->pc     = pc;
47674ae819SStefano Zampini   *fetidpmat_ctx = newctx;
48674ae819SStefano Zampini   PetscFunctionReturn(0);
49674ae819SStefano Zampini }
50674ae819SStefano Zampini 
519371c9d4SSatish Balay PetscErrorCode PCBDDCCreateFETIDPPCContext(PC pc, FETIDPPC_ctx *fetidppc_ctx) {
52674ae819SStefano Zampini   FETIDPPC_ctx newctx;
53674ae819SStefano Zampini 
54674ae819SStefano Zampini   PetscFunctionBegin;
559566063dSJacob Faibussowitsch   PetscCall(PetscNew(&newctx));
56674ae819SStefano Zampini   /* increase the reference count for BDDC preconditioner */
579566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)pc));
58674ae819SStefano Zampini   newctx->pc    = pc;
59674ae819SStefano Zampini   *fetidppc_ctx = newctx;
60674ae819SStefano Zampini   PetscFunctionReturn(0);
61674ae819SStefano Zampini }
62674ae819SStefano Zampini 
639371c9d4SSatish Balay PetscErrorCode PCBDDCDestroyFETIDPMat(Mat A) {
64674ae819SStefano Zampini   FETIDPMat_ctx mat_ctx;
65674ae819SStefano Zampini 
66674ae819SStefano Zampini   PetscFunctionBegin;
679566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(A, &mat_ctx));
689566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&mat_ctx->lambda_local));
699566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&mat_ctx->temp_solution_D));
709566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&mat_ctx->temp_solution_B));
719566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&mat_ctx->B_delta));
729566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&mat_ctx->B_Ddelta));
739566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&mat_ctx->B_BB));
749566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&mat_ctx->B_BI));
759566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&mat_ctx->Bt_BB));
769566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&mat_ctx->Bt_BI));
779566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&mat_ctx->C));
789566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&mat_ctx->rhs_flip));
799566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&mat_ctx->vP));
809566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&mat_ctx->xPg));
819566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&mat_ctx->yPg));
829566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&mat_ctx->l2g_lambda));
839566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&mat_ctx->l2g_lambda_only));
849566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&mat_ctx->l2g_p));
859566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&mat_ctx->g2g_p));
869566063dSJacob Faibussowitsch   PetscCall(PCDestroy(&mat_ctx->pc)); /* decrease PCBDDC reference count */
879566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&mat_ctx->pressure));
889566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&mat_ctx->lagrange));
899566063dSJacob Faibussowitsch   PetscCall(PetscFree(mat_ctx));
90674ae819SStefano Zampini   PetscFunctionReturn(0);
91674ae819SStefano Zampini }
92674ae819SStefano Zampini 
939371c9d4SSatish Balay PetscErrorCode PCBDDCDestroyFETIDPPC(PC pc) {
94674ae819SStefano Zampini   FETIDPPC_ctx pc_ctx;
95674ae819SStefano Zampini 
96674ae819SStefano Zampini   PetscFunctionBegin;
979566063dSJacob Faibussowitsch   PetscCall(PCShellGetContext(pc, &pc_ctx));
989566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pc_ctx->lambda_local));
999566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pc_ctx->B_Ddelta));
1009566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&pc_ctx->l2g_lambda));
1019566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pc_ctx->S_j));
1029566063dSJacob Faibussowitsch   PetscCall(PCDestroy(&pc_ctx->pc)); /* decrease PCBDDC reference count */
1039566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pc_ctx->xPg));
1049566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pc_ctx->yPg));
1059566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc_ctx));
106674ae819SStefano Zampini   PetscFunctionReturn(0);
107674ae819SStefano Zampini }
108674ae819SStefano Zampini 
1099371c9d4SSatish Balay PetscErrorCode PCBDDCSetupFETIDPMatContext(FETIDPMat_ctx fetidpmat_ctx) {
110674ae819SStefano Zampini   PC_IS          *pcis      = (PC_IS *)fetidpmat_ctx->pc->data;
111674ae819SStefano Zampini   PC_BDDC        *pcbddc    = (PC_BDDC *)fetidpmat_ctx->pc->data;
112674ae819SStefano Zampini   PCBDDCGraph     mat_graph = pcbddc->mat_graph;
113674ae819SStefano Zampini   Mat_IS         *matis     = (Mat_IS *)fetidpmat_ctx->pc->pmat->data;
114674ae819SStefano Zampini   MPI_Comm        comm;
11525d06dbeSstefano_zampini   Mat             ScalingMat, BD1, BD2;
116457d4a33Sstefano_zampini   Vec             fetidp_global;
117674ae819SStefano Zampini   IS              IS_l2g_lambda;
118a1c0d0daSstefano_zampini   IS              subset, subset_mult, subset_n, isvert;
119674ae819SStefano Zampini   PetscBool       skip_node, fully_redundant;
120674ae819SStefano Zampini   PetscInt        i, j, k, s, n_boundary_dofs, n_global_lambda, n_vertices, partial_sum;
121dc456d91SStefano Zampini   PetscInt        cum, n_local_lambda, n_lambda_for_dof, dual_size, n_neg_values, n_pos_values;
12276ec1555SBarry Smith   PetscMPIInt     rank, size, buf_size, neigh;
123674ae819SStefano Zampini   PetscScalar     scalar_value;
124a1c0d0daSstefano_zampini   const PetscInt *vertex_indices;
125dc456d91SStefano Zampini   PetscInt       *dual_dofs_boundary_indices, *aux_local_numbering_1;
126dc456d91SStefano Zampini   const PetscInt *aux_global_numbering;
127674ae819SStefano Zampini   PetscInt       *aux_sums, *cols_B_delta, *l2g_indices;
128674ae819SStefano Zampini   PetscScalar    *array, *scaling_factors, *vals_B_delta;
12925d06dbeSstefano_zampini   PetscScalar   **all_factors;
130674ae819SStefano Zampini   PetscInt       *aux_local_numbering_2;
131457d4a33Sstefano_zampini   PetscLayout     llay;
132a1c0d0daSstefano_zampini 
133457d4a33Sstefano_zampini   /* saddlepoint */
134457d4a33Sstefano_zampini   ISLocalToGlobalMapping l2gmap_p;
135457d4a33Sstefano_zampini   PetscLayout            play;
136457d4a33Sstefano_zampini   IS                     gP, pP;
137457d4a33Sstefano_zampini   PetscInt               nPl, nPg, nPgl;
138674ae819SStefano Zampini 
139674ae819SStefano Zampini   PetscFunctionBegin;
1409566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)(fetidpmat_ctx->pc), &comm));
1419566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
1429566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
143674ae819SStefano Zampini 
144457d4a33Sstefano_zampini   /* saddlepoint */
145457d4a33Sstefano_zampini   nPl      = 0;
146457d4a33Sstefano_zampini   nPg      = 0;
147457d4a33Sstefano_zampini   nPgl     = 0;
148457d4a33Sstefano_zampini   gP       = NULL;
149457d4a33Sstefano_zampini   pP       = NULL;
150457d4a33Sstefano_zampini   l2gmap_p = NULL;
151457d4a33Sstefano_zampini   play     = NULL;
1529566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)fetidpmat_ctx->pc, "__KSPFETIDP_pP", (PetscObject *)&pP));
153022d8d2bSstefano_zampini   if (pP) { /* saddle point */
154457d4a33Sstefano_zampini     /* subdomain pressures in global numbering */
1559566063dSJacob Faibussowitsch     PetscCall(PetscObjectQuery((PetscObject)fetidpmat_ctx->pc, "__KSPFETIDP_gP", (PetscObject *)&gP));
15628b400f6SJacob Faibussowitsch     PetscCheck(gP, PETSC_COMM_SELF, PETSC_ERR_PLIB, "gP not present");
1579566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(gP, &nPl));
1589566063dSJacob Faibussowitsch     PetscCall(VecCreate(PETSC_COMM_SELF, &fetidpmat_ctx->vP));
1599566063dSJacob Faibussowitsch     PetscCall(VecSetSizes(fetidpmat_ctx->vP, nPl, nPl));
1609566063dSJacob Faibussowitsch     PetscCall(VecSetType(fetidpmat_ctx->vP, VECSTANDARD));
1619566063dSJacob Faibussowitsch     PetscCall(VecSetUp(fetidpmat_ctx->vP));
162457d4a33Sstefano_zampini 
163e1214c54Sstefano_zampini     /* pressure matrix */
1649566063dSJacob Faibussowitsch     PetscCall(PetscObjectQuery((PetscObject)fetidpmat_ctx->pc, "__KSPFETIDP_C", (PetscObject *)&fetidpmat_ctx->C));
165e1214c54Sstefano_zampini     if (!fetidpmat_ctx->C) { /* null pressure block, compute layout and global numbering for pressures */
166457d4a33Sstefano_zampini       IS Pg;
167457d4a33Sstefano_zampini 
1689566063dSJacob Faibussowitsch       PetscCall(ISRenumber(gP, NULL, &nPg, &Pg));
1699566063dSJacob Faibussowitsch       PetscCall(ISLocalToGlobalMappingCreateIS(Pg, &l2gmap_p));
1709566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&Pg));
1719566063dSJacob Faibussowitsch       PetscCall(PetscLayoutCreate(comm, &play));
1729566063dSJacob Faibussowitsch       PetscCall(PetscLayoutSetBlockSize(play, 1));
1739566063dSJacob Faibussowitsch       PetscCall(PetscLayoutSetSize(play, nPg));
1749566063dSJacob Faibussowitsch       PetscCall(ISGetLocalSize(pP, &nPgl));
1759566063dSJacob Faibussowitsch       PetscCall(PetscLayoutSetLocalSize(play, nPgl));
1769566063dSJacob Faibussowitsch       PetscCall(PetscLayoutSetUp(play));
177457d4a33Sstefano_zampini     } else {
1789566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)fetidpmat_ctx->C));
1799566063dSJacob Faibussowitsch       PetscCall(MatISGetLocalToGlobalMapping(fetidpmat_ctx->C, &l2gmap_p, NULL));
1809566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)l2gmap_p));
1819566063dSJacob Faibussowitsch       PetscCall(MatGetSize(fetidpmat_ctx->C, &nPg, NULL));
1829566063dSJacob Faibussowitsch       PetscCall(MatGetLocalSize(fetidpmat_ctx->C, NULL, &nPgl));
1839566063dSJacob Faibussowitsch       PetscCall(MatGetLayouts(fetidpmat_ctx->C, NULL, &llay));
1849566063dSJacob Faibussowitsch       PetscCall(PetscLayoutReference(llay, &play));
185457d4a33Sstefano_zampini     }
1869566063dSJacob Faibussowitsch     PetscCall(VecCreateMPIWithArray(comm, 1, nPgl, nPg, NULL, &fetidpmat_ctx->xPg));
1879566063dSJacob Faibussowitsch     PetscCall(VecCreateMPIWithArray(comm, 1, nPgl, nPg, NULL, &fetidpmat_ctx->yPg));
188457d4a33Sstefano_zampini 
189e1214c54Sstefano_zampini     /* import matrices for pressures coupling */
1909566063dSJacob Faibussowitsch     PetscCall(PetscObjectQuery((PetscObject)fetidpmat_ctx->pc, "__KSPFETIDP_B_BI", (PetscObject *)&fetidpmat_ctx->B_BI));
19128b400f6SJacob Faibussowitsch     PetscCheck(fetidpmat_ctx->B_BI, PETSC_COMM_SELF, PETSC_ERR_PLIB, "B_BI not present");
1929566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)fetidpmat_ctx->B_BI));
193457d4a33Sstefano_zampini 
1949566063dSJacob Faibussowitsch     PetscCall(PetscObjectQuery((PetscObject)fetidpmat_ctx->pc, "__KSPFETIDP_B_BB", (PetscObject *)&fetidpmat_ctx->B_BB));
19528b400f6SJacob Faibussowitsch     PetscCheck(fetidpmat_ctx->B_BB, PETSC_COMM_SELF, PETSC_ERR_PLIB, "B_BB not present");
1969566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)fetidpmat_ctx->B_BB));
197457d4a33Sstefano_zampini 
1989566063dSJacob Faibussowitsch     PetscCall(PetscObjectQuery((PetscObject)fetidpmat_ctx->pc, "__KSPFETIDP_Bt_BI", (PetscObject *)&fetidpmat_ctx->Bt_BI));
19928b400f6SJacob Faibussowitsch     PetscCheck(fetidpmat_ctx->Bt_BI, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Bt_BI not present");
2009566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)fetidpmat_ctx->Bt_BI));
201457d4a33Sstefano_zampini 
2029566063dSJacob Faibussowitsch     PetscCall(PetscObjectQuery((PetscObject)fetidpmat_ctx->pc, "__KSPFETIDP_Bt_BB", (PetscObject *)&fetidpmat_ctx->Bt_BB));
20328b400f6SJacob Faibussowitsch     PetscCheck(fetidpmat_ctx->Bt_BB, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Bt_BB not present");
2049566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)fetidpmat_ctx->Bt_BB));
2056cc1294bSstefano_zampini 
2069566063dSJacob Faibussowitsch     PetscCall(PetscObjectQuery((PetscObject)fetidpmat_ctx->pc, "__KSPFETIDP_flip", (PetscObject *)&fetidpmat_ctx->rhs_flip));
2071baa6e33SBarry Smith     if (fetidpmat_ctx->rhs_flip) PetscCall(PetscObjectReference((PetscObject)fetidpmat_ctx->rhs_flip));
208457d4a33Sstefano_zampini   }
209457d4a33Sstefano_zampini 
210674ae819SStefano Zampini   /* Default type of lagrange multipliers is non-redundant */
211329cd39dSStefano Zampini   fully_redundant = fetidpmat_ctx->fully_redundant;
212674ae819SStefano Zampini 
213674ae819SStefano Zampini   /* Evaluate local and global number of lagrange multipliers */
2149566063dSJacob Faibussowitsch   PetscCall(VecSet(pcis->vec1_N, 0.0));
215674ae819SStefano Zampini   n_local_lambda  = 0;
216674ae819SStefano Zampini   partial_sum     = 0;
217674ae819SStefano Zampini   n_boundary_dofs = 0;
218674ae819SStefano Zampini   s               = 0;
219a1c0d0daSstefano_zampini 
220674ae819SStefano Zampini   /* Get Vertices used to define the BDDC */
2219566063dSJacob Faibussowitsch   PetscCall(PCBDDCGraphGetCandidatesIS(pcbddc->mat_graph, NULL, NULL, NULL, NULL, &isvert));
2229566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(isvert, &n_vertices));
2239566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(isvert, &vertex_indices));
224a1c0d0daSstefano_zampini 
225674ae819SStefano Zampini   dual_size = pcis->n_B - n_vertices;
2269566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(dual_size, &dual_dofs_boundary_indices));
2279566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(dual_size, &aux_local_numbering_1));
2289566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(dual_size, &aux_local_numbering_2));
229674ae819SStefano Zampini 
2309566063dSJacob Faibussowitsch   PetscCall(VecGetArray(pcis->vec1_N, &array));
231674ae819SStefano Zampini   for (i = 0; i < pcis->n; i++) {
232674ae819SStefano Zampini     j = mat_graph->count[i]; /* RECALL: mat_graph->count[i] does not count myself */
2332d456bbaSstefano_zampini     if (j > 0) n_boundary_dofs++;
234674ae819SStefano Zampini     skip_node = PETSC_FALSE;
235674ae819SStefano Zampini     if (s < n_vertices && vertex_indices[s] == i) { /* it works for a sorted set of vertices */
236674ae819SStefano Zampini       skip_node = PETSC_TRUE;
237674ae819SStefano Zampini       s++;
238674ae819SStefano Zampini     }
2392d456bbaSstefano_zampini     if (j < 1) skip_node = PETSC_TRUE;
2402d456bbaSstefano_zampini     if (mat_graph->special_dof[i] == PCBDDCGRAPH_DIRICHLET_MARK) skip_node = PETSC_TRUE;
241674ae819SStefano Zampini     if (!skip_node) {
242674ae819SStefano Zampini       if (fully_redundant) {
243674ae819SStefano Zampini         /* fully redundant set of lagrange multipliers */
244674ae819SStefano Zampini         n_lambda_for_dof = (j * (j + 1)) / 2;
245674ae819SStefano Zampini       } else {
246674ae819SStefano Zampini         n_lambda_for_dof = j;
247674ae819SStefano Zampini       }
248674ae819SStefano Zampini       n_local_lambda += j;
249674ae819SStefano Zampini       /* needed to evaluate global number of lagrange multipliers */
250674ae819SStefano Zampini       array[i]                                = (1.0 * n_lambda_for_dof) / (j + 1.0); /* already scaled for the next global sum */
251674ae819SStefano Zampini       /* store some data needed */
252674ae819SStefano Zampini       dual_dofs_boundary_indices[partial_sum] = n_boundary_dofs - 1;
253674ae819SStefano Zampini       aux_local_numbering_1[partial_sum]      = i;
254674ae819SStefano Zampini       aux_local_numbering_2[partial_sum]      = n_lambda_for_dof;
255674ae819SStefano Zampini       partial_sum++;
256674ae819SStefano Zampini     }
257674ae819SStefano Zampini   }
2589566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(pcis->vec1_N, &array));
2599566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isvert, &vertex_indices));
2609566063dSJacob Faibussowitsch   PetscCall(PCBDDCGraphRestoreCandidatesIS(pcbddc->mat_graph, NULL, NULL, NULL, NULL, &isvert));
2612d456bbaSstefano_zampini   dual_size = partial_sum;
262674ae819SStefano Zampini 
263674ae819SStefano Zampini   /* compute global ordering of lagrange multipliers and associate l2g map */
2649566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(comm, partial_sum, aux_local_numbering_1, PETSC_COPY_VALUES, &subset_n));
2659566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingApplyIS(pcis->mapping, subset_n, &subset));
2669566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&subset_n));
2679566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(comm, partial_sum, aux_local_numbering_2, PETSC_OWN_POINTER, &subset_mult));
2689566063dSJacob Faibussowitsch   PetscCall(ISRenumber(subset, subset_mult, &fetidpmat_ctx->n_lambda, &subset_n));
2699566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&subset));
2703d996552SStefano Zampini 
27176bd3646SJed Brown   if (PetscDefined(USE_DEBUG)) {
2729566063dSJacob Faibussowitsch     PetscCall(VecSet(pcis->vec1_global, 0.0));
2739566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(matis->rctx, pcis->vec1_N, pcis->vec1_global, ADD_VALUES, SCATTER_REVERSE));
2749566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(matis->rctx, pcis->vec1_N, pcis->vec1_global, ADD_VALUES, SCATTER_REVERSE));
2759566063dSJacob Faibussowitsch     PetscCall(VecSum(pcis->vec1_global, &scalar_value));
2764fe826edSStefano Zampini     i = (PetscInt)PetscRealPart(scalar_value);
27763a3b9bcSJacob 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);
27876bd3646SJed Brown   }
279674ae819SStefano Zampini 
280674ae819SStefano Zampini   /* init data for scaling factors exchange */
28125d06dbeSstefano_zampini   if (!pcbddc->use_deluxe_scaling) {
28225d06dbeSstefano_zampini     PetscInt    *ptrs_buffer, neigh_position;
28325d06dbeSstefano_zampini     PetscScalar *send_buffer, *recv_buffer;
28425d06dbeSstefano_zampini     MPI_Request *send_reqs, *recv_reqs;
28525d06dbeSstefano_zampini 
286674ae819SStefano Zampini     partial_sum = 0;
2879566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(pcis->n_neigh, &ptrs_buffer));
2889566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(PetscMax(pcis->n_neigh - 1, 0), &send_reqs));
2899566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(PetscMax(pcis->n_neigh - 1, 0), &recv_reqs));
2909566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(pcis->n + 1, &all_factors));
2914b2aedd3SStefano Zampini     if (pcis->n_neigh > 0) ptrs_buffer[0] = 0;
292674ae819SStefano Zampini     for (i = 1; i < pcis->n_neigh; i++) {
293674ae819SStefano Zampini       partial_sum += pcis->n_shared[i];
294674ae819SStefano Zampini       ptrs_buffer[i] = ptrs_buffer[i - 1] + pcis->n_shared[i];
295674ae819SStefano Zampini     }
2969566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(partial_sum, &send_buffer));
2979566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(partial_sum, &recv_buffer));
2989566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(partial_sum, &all_factors[0]));
299674ae819SStefano Zampini     for (i = 0; i < pcis->n - 1; i++) {
300674ae819SStefano Zampini       j                  = mat_graph->count[i];
301674ae819SStefano Zampini       all_factors[i + 1] = all_factors[i] + j;
302674ae819SStefano Zampini     }
30325d06dbeSstefano_zampini 
304674ae819SStefano Zampini     /* scatter B scaling to N vec */
3059566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(pcis->N_to_B, pcis->D, pcis->vec1_N, INSERT_VALUES, SCATTER_REVERSE));
3069566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(pcis->N_to_B, pcis->D, pcis->vec1_N, INSERT_VALUES, SCATTER_REVERSE));
307674ae819SStefano Zampini     /* communications */
3089566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(pcis->vec1_N, (const PetscScalar **)&array));
309674ae819SStefano Zampini     for (i = 1; i < pcis->n_neigh; i++) {
3109371c9d4SSatish Balay       for (j = 0; j < pcis->n_shared[i]; j++) { send_buffer[ptrs_buffer[i - 1] + j] = array[pcis->shared[i][j]]; }
3119566063dSJacob Faibussowitsch       PetscCall(PetscMPIIntCast(ptrs_buffer[i] - ptrs_buffer[i - 1], &buf_size));
3129566063dSJacob Faibussowitsch       PetscCall(PetscMPIIntCast(pcis->neigh[i], &neigh));
3139566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Isend(&send_buffer[ptrs_buffer[i - 1]], buf_size, MPIU_SCALAR, neigh, 0, comm, &send_reqs[i - 1]));
3149566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Irecv(&recv_buffer[ptrs_buffer[i - 1]], buf_size, MPIU_SCALAR, neigh, 0, comm, &recv_reqs[i - 1]));
315674ae819SStefano Zampini     }
3169566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(pcis->vec1_N, (const PetscScalar **)&array));
317*48a46eb9SPierre Jolivet     if (pcis->n_neigh > 0) PetscCallMPI(MPI_Waitall(pcis->n_neigh - 1, recv_reqs, MPI_STATUSES_IGNORE));
318674ae819SStefano Zampini     /* put values in correct places */
319674ae819SStefano Zampini     for (i = 1; i < pcis->n_neigh; i++) {
320674ae819SStefano Zampini       for (j = 0; j < pcis->n_shared[i]; j++) {
321674ae819SStefano Zampini         k              = pcis->shared[i][j];
322674ae819SStefano Zampini         neigh_position = 0;
323674ae819SStefano Zampini         while (mat_graph->neighbours_set[k][neigh_position] != pcis->neigh[i]) { neigh_position++; }
324674ae819SStefano Zampini         all_factors[k][neigh_position] = recv_buffer[ptrs_buffer[i - 1] + j];
325674ae819SStefano Zampini       }
326674ae819SStefano Zampini     }
327*48a46eb9SPierre Jolivet     if (pcis->n_neigh > 0) PetscCallMPI(MPI_Waitall(pcis->n_neigh - 1, send_reqs, MPI_STATUSES_IGNORE));
3289566063dSJacob Faibussowitsch     PetscCall(PetscFree(send_reqs));
3299566063dSJacob Faibussowitsch     PetscCall(PetscFree(recv_reqs));
3309566063dSJacob Faibussowitsch     PetscCall(PetscFree(send_buffer));
3319566063dSJacob Faibussowitsch     PetscCall(PetscFree(recv_buffer));
3329566063dSJacob Faibussowitsch     PetscCall(PetscFree(ptrs_buffer));
33325d06dbeSstefano_zampini   }
334674ae819SStefano Zampini 
335674ae819SStefano Zampini   /* Compute B and B_delta (local actions) */
3369566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(pcis->n_neigh, &aux_sums));
3379566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n_local_lambda, &l2g_indices));
3389566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n_local_lambda, &vals_B_delta));
3399566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n_local_lambda, &cols_B_delta));
34025d06dbeSstefano_zampini   if (!pcbddc->use_deluxe_scaling) {
3419566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n_local_lambda, &scaling_factors));
34225d06dbeSstefano_zampini   } else {
34325d06dbeSstefano_zampini     scaling_factors = NULL;
34425d06dbeSstefano_zampini     all_factors     = NULL;
34525d06dbeSstefano_zampini   }
3469566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(subset_n, &aux_global_numbering));
347674ae819SStefano Zampini   partial_sum = 0;
348dc456d91SStefano Zampini   cum         = 0;
349674ae819SStefano Zampini   for (i = 0; i < dual_size; i++) {
350dc456d91SStefano Zampini     n_global_lambda = aux_global_numbering[cum];
351674ae819SStefano Zampini     j               = mat_graph->count[aux_local_numbering_1[i]];
352674ae819SStefano Zampini     aux_sums[0]     = 0;
3539371c9d4SSatish Balay     for (s = 1; s < j; s++) { aux_sums[s] = aux_sums[s - 1] + j - s + 1; }
35425d06dbeSstefano_zampini     if (all_factors) array = all_factors[aux_local_numbering_1[i]];
355674ae819SStefano Zampini     n_neg_values = 0;
3562a7da448SStefano Zampini     while (n_neg_values < j && mat_graph->neighbours_set[aux_local_numbering_1[i]][n_neg_values] < rank) { n_neg_values++; }
357674ae819SStefano Zampini     n_pos_values = j - n_neg_values;
358674ae819SStefano Zampini     if (fully_redundant) {
359674ae819SStefano Zampini       for (s = 0; s < n_neg_values; s++) {
360674ae819SStefano Zampini         l2g_indices[partial_sum + s]  = aux_sums[s] + n_neg_values - s - 1 + n_global_lambda;
361674ae819SStefano Zampini         cols_B_delta[partial_sum + s] = dual_dofs_boundary_indices[i];
362674ae819SStefano Zampini         vals_B_delta[partial_sum + s] = -1.0;
36325d06dbeSstefano_zampini         if (!pcbddc->use_deluxe_scaling) scaling_factors[partial_sum + s] = array[s];
364674ae819SStefano Zampini       }
365674ae819SStefano Zampini       for (s = 0; s < n_pos_values; s++) {
366674ae819SStefano Zampini         l2g_indices[partial_sum + s + n_neg_values]  = aux_sums[n_neg_values] + s + n_global_lambda;
367674ae819SStefano Zampini         cols_B_delta[partial_sum + s + n_neg_values] = dual_dofs_boundary_indices[i];
368674ae819SStefano Zampini         vals_B_delta[partial_sum + s + n_neg_values] = 1.0;
36925d06dbeSstefano_zampini         if (!pcbddc->use_deluxe_scaling) scaling_factors[partial_sum + s + n_neg_values] = array[s + n_neg_values];
370674ae819SStefano Zampini       }
371674ae819SStefano Zampini       partial_sum += j;
372674ae819SStefano Zampini     } else {
373674ae819SStefano Zampini       /* l2g_indices and default cols and vals of B_delta */
374674ae819SStefano Zampini       for (s = 0; s < j; s++) {
375674ae819SStefano Zampini         l2g_indices[partial_sum + s]  = n_global_lambda + s;
376674ae819SStefano Zampini         cols_B_delta[partial_sum + s] = dual_dofs_boundary_indices[i];
377674ae819SStefano Zampini         vals_B_delta[partial_sum + s] = 0.0;
378674ae819SStefano Zampini       }
379674ae819SStefano Zampini       /* B_delta */
380674ae819SStefano Zampini       if (n_neg_values > 0) { /* there's a rank next to me to the left */
381674ae819SStefano Zampini         vals_B_delta[partial_sum + n_neg_values - 1] = -1.0;
382674ae819SStefano Zampini       }
383674ae819SStefano Zampini       if (n_neg_values < j) { /* there's a rank next to me to the right */
384674ae819SStefano Zampini         vals_B_delta[partial_sum + n_neg_values] = 1.0;
385674ae819SStefano Zampini       }
386674ae819SStefano Zampini       /* scaling as in Klawonn-Widlund 1999 */
38725d06dbeSstefano_zampini       if (!pcbddc->use_deluxe_scaling) {
388674ae819SStefano Zampini         for (s = 0; s < n_neg_values; s++) {
389674ae819SStefano Zampini           scalar_value = 0.0;
39025d06dbeSstefano_zampini           for (k = 0; k < s + 1; k++) scalar_value += array[k];
391674ae819SStefano Zampini           scaling_factors[partial_sum + s] = -scalar_value;
392674ae819SStefano Zampini         }
393674ae819SStefano Zampini         for (s = 0; s < n_pos_values; s++) {
394674ae819SStefano Zampini           scalar_value = 0.0;
39525d06dbeSstefano_zampini           for (k = s + n_neg_values; k < j; k++) scalar_value += array[k];
396674ae819SStefano Zampini           scaling_factors[partial_sum + s + n_neg_values] = scalar_value;
397674ae819SStefano Zampini         }
39825d06dbeSstefano_zampini       }
399674ae819SStefano Zampini       partial_sum += j;
400674ae819SStefano Zampini     }
401dc456d91SStefano Zampini     cum += aux_local_numbering_2[i];
402674ae819SStefano Zampini   }
4039566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(subset_n, &aux_global_numbering));
4049566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&subset_mult));
4059566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&subset_n));
4069566063dSJacob Faibussowitsch   PetscCall(PetscFree(aux_sums));
4079566063dSJacob Faibussowitsch   PetscCall(PetscFree(aux_local_numbering_1));
4089566063dSJacob Faibussowitsch   PetscCall(PetscFree(dual_dofs_boundary_indices));
40925d06dbeSstefano_zampini   if (all_factors) {
4109566063dSJacob Faibussowitsch     PetscCall(PetscFree(all_factors[0]));
4119566063dSJacob Faibussowitsch     PetscCall(PetscFree(all_factors));
41225d06dbeSstefano_zampini   }
413674ae819SStefano Zampini 
414674ae819SStefano Zampini   /* Create local part of B_delta */
4159566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_SELF, &fetidpmat_ctx->B_delta));
4169566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(fetidpmat_ctx->B_delta, n_local_lambda, pcis->n_B, n_local_lambda, pcis->n_B));
4179566063dSJacob Faibussowitsch   PetscCall(MatSetType(fetidpmat_ctx->B_delta, MATSEQAIJ));
4189566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(fetidpmat_ctx->B_delta, 1, NULL));
4199566063dSJacob Faibussowitsch   PetscCall(MatSetOption(fetidpmat_ctx->B_delta, MAT_IGNORE_ZERO_ENTRIES, PETSC_TRUE));
420*48a46eb9SPierre 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));
4219566063dSJacob Faibussowitsch   PetscCall(PetscFree(vals_B_delta));
4229566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(fetidpmat_ctx->B_delta, MAT_FINAL_ASSEMBLY));
4239566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(fetidpmat_ctx->B_delta, MAT_FINAL_ASSEMBLY));
424674ae819SStefano Zampini 
42525d06dbeSstefano_zampini   BD1 = NULL;
42625d06dbeSstefano_zampini   BD2 = NULL;
427674ae819SStefano Zampini   if (fully_redundant) {
42828b400f6SJacob Faibussowitsch     PetscCheck(!pcbddc->use_deluxe_scaling, comm, PETSC_ERR_SUP, "Deluxe FETIDP with fully-redundant multipliers to be implemented");
4299566063dSJacob Faibussowitsch     PetscCall(MatCreate(PETSC_COMM_SELF, &ScalingMat));
4309566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(ScalingMat, n_local_lambda, n_local_lambda, n_local_lambda, n_local_lambda));
4319566063dSJacob Faibussowitsch     PetscCall(MatSetType(ScalingMat, MATSEQAIJ));
4329566063dSJacob Faibussowitsch     PetscCall(MatSeqAIJSetPreallocation(ScalingMat, 1, NULL));
433*48a46eb9SPierre Jolivet     for (i = 0; i < n_local_lambda; i++) PetscCall(MatSetValue(ScalingMat, i, i, scaling_factors[i], INSERT_VALUES));
4349566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(ScalingMat, MAT_FINAL_ASSEMBLY));
4359566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(ScalingMat, MAT_FINAL_ASSEMBLY));
4369566063dSJacob Faibussowitsch     PetscCall(MatMatMult(ScalingMat, fetidpmat_ctx->B_delta, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &fetidpmat_ctx->B_Ddelta));
4379566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&ScalingMat));
438674ae819SStefano Zampini   } else {
4399566063dSJacob Faibussowitsch     PetscCall(MatCreate(PETSC_COMM_SELF, &fetidpmat_ctx->B_Ddelta));
4409566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(fetidpmat_ctx->B_Ddelta, n_local_lambda, pcis->n_B, n_local_lambda, pcis->n_B));
441524221abSStefano Zampini     if (!pcbddc->use_deluxe_scaling || !pcbddc->sub_schurs) {
4429566063dSJacob Faibussowitsch       PetscCall(MatSetType(fetidpmat_ctx->B_Ddelta, MATSEQAIJ));
4439566063dSJacob Faibussowitsch       PetscCall(MatSeqAIJSetPreallocation(fetidpmat_ctx->B_Ddelta, 1, NULL));
444*48a46eb9SPierre 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));
4459566063dSJacob Faibussowitsch       PetscCall(MatAssemblyBegin(fetidpmat_ctx->B_Ddelta, MAT_FINAL_ASSEMBLY));
4469566063dSJacob Faibussowitsch       PetscCall(MatAssemblyEnd(fetidpmat_ctx->B_Ddelta, MAT_FINAL_ASSEMBLY));
44725d06dbeSstefano_zampini     } else {
44825d06dbeSstefano_zampini       /* scaling as in Klawonn-Widlund 1999 */
44925d06dbeSstefano_zampini       PCBDDCDeluxeScaling deluxe_ctx = pcbddc->deluxe_ctx;
45025d06dbeSstefano_zampini       PCBDDCSubSchurs     sub_schurs = pcbddc->sub_schurs;
45125d06dbeSstefano_zampini       Mat                 T;
45225d06dbeSstefano_zampini       PetscScalar        *W, lwork, *Bwork;
453451a39c7SStefano Zampini       const PetscInt     *idxs = NULL;
45425d06dbeSstefano_zampini       PetscInt            cum, mss, *nnz;
45525d06dbeSstefano_zampini       PetscBLASInt       *pivots, B_lwork, B_N, B_ierr;
45625d06dbeSstefano_zampini 
45728b400f6SJacob Faibussowitsch       PetscCheck(pcbddc->deluxe_singlemat, comm, PETSC_ERR_USER, "Cannot compute B_Ddelta! rerun with -pc_bddc_deluxe_singlemat");
45825d06dbeSstefano_zampini       mss = 0;
4599566063dSJacob Faibussowitsch       PetscCall(PetscCalloc1(pcis->n_B, &nnz));
46025d06dbeSstefano_zampini       if (sub_schurs->is_Ej_all) {
4619566063dSJacob Faibussowitsch         PetscCall(ISGetIndices(sub_schurs->is_Ej_all, &idxs));
46225d06dbeSstefano_zampini         for (i = 0, cum = 0; i < sub_schurs->n_subs; i++) {
46325d06dbeSstefano_zampini           PetscInt subset_size;
46425d06dbeSstefano_zampini 
4659566063dSJacob Faibussowitsch           PetscCall(ISGetLocalSize(sub_schurs->is_subs[i], &subset_size));
46625d06dbeSstefano_zampini           for (j = 0; j < subset_size; j++) nnz[idxs[j + cum]] = subset_size;
46725d06dbeSstefano_zampini           mss = PetscMax(mss, subset_size);
46825d06dbeSstefano_zampini           cum += subset_size;
46925d06dbeSstefano_zampini         }
47025d06dbeSstefano_zampini       }
4719566063dSJacob Faibussowitsch       PetscCall(MatCreate(PETSC_COMM_SELF, &T));
4729566063dSJacob Faibussowitsch       PetscCall(MatSetSizes(T, pcis->n_B, pcis->n_B, pcis->n_B, pcis->n_B));
4739566063dSJacob Faibussowitsch       PetscCall(MatSetType(T, MATSEQAIJ));
4749566063dSJacob Faibussowitsch       PetscCall(MatSeqAIJSetPreallocation(T, 0, nnz));
4759566063dSJacob Faibussowitsch       PetscCall(PetscFree(nnz));
47625d06dbeSstefano_zampini 
47725d06dbeSstefano_zampini       /* workspace allocation */
4787f00194aSstefano_zampini       B_lwork = 0;
4797f00194aSstefano_zampini       if (mss) {
480be8bcbd0Sstefano_zampini         PetscScalar dummy = 1;
481be8bcbd0Sstefano_zampini 
48225d06dbeSstefano_zampini         B_lwork = -1;
4839566063dSJacob Faibussowitsch         PetscCall(PetscBLASIntCast(mss, &B_N));
4849566063dSJacob Faibussowitsch         PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
485792fecdfSBarry Smith         PetscCallBLAS("LAPACKgetri", LAPACKgetri_(&B_N, &dummy, &B_N, &B_N, &lwork, &B_lwork, &B_ierr));
4869566063dSJacob Faibussowitsch         PetscCall(PetscFPTrapPop());
48728b400f6SJacob Faibussowitsch         PetscCheck(!B_ierr, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in query to GETRI Lapack routine %d", (int)B_ierr);
4889566063dSJacob Faibussowitsch         PetscCall(PetscBLASIntCast((PetscInt)PetscRealPart(lwork), &B_lwork));
4897f00194aSstefano_zampini       }
4909566063dSJacob Faibussowitsch       PetscCall(PetscMalloc3(mss * mss, &W, mss, &pivots, B_lwork, &Bwork));
49125d06dbeSstefano_zampini 
49225d06dbeSstefano_zampini       for (i = 0, cum = 0; i < sub_schurs->n_subs; i++) {
4931683a169SBarry Smith         const PetscScalar *M;
49425d06dbeSstefano_zampini         PetscInt           subset_size;
49525d06dbeSstefano_zampini 
4969566063dSJacob Faibussowitsch         PetscCall(ISGetLocalSize(sub_schurs->is_subs[i], &subset_size));
4979566063dSJacob Faibussowitsch         PetscCall(PetscBLASIntCast(subset_size, &B_N));
4989566063dSJacob Faibussowitsch         PetscCall(MatDenseGetArrayRead(deluxe_ctx->seq_mat[i], &M));
4999566063dSJacob Faibussowitsch         PetscCall(PetscArraycpy(W, M, subset_size * subset_size));
5009566063dSJacob Faibussowitsch         PetscCall(MatDenseRestoreArrayRead(deluxe_ctx->seq_mat[i], &M));
5019566063dSJacob Faibussowitsch         PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
502792fecdfSBarry Smith         PetscCallBLAS("LAPACKgetrf", LAPACKgetrf_(&B_N, &B_N, W, &B_N, pivots, &B_ierr));
50328b400f6SJacob Faibussowitsch         PetscCheck(!B_ierr, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in GETRF Lapack routine %d", (int)B_ierr);
504792fecdfSBarry Smith         PetscCallBLAS("LAPACKgetri", LAPACKgetri_(&B_N, W, &B_N, pivots, Bwork, &B_lwork, &B_ierr));
50528b400f6SJacob Faibussowitsch         PetscCheck(!B_ierr, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in GETRI Lapack routine %d", (int)B_ierr);
5069566063dSJacob Faibussowitsch         PetscCall(PetscFPTrapPop());
507451a39c7SStefano Zampini         /* silent static analyzer */
50828b400f6SJacob Faibussowitsch         PetscCheck(idxs, PETSC_COMM_SELF, PETSC_ERR_PLIB, "IDXS not present");
5099566063dSJacob Faibussowitsch         PetscCall(MatSetValues(T, subset_size, idxs + cum, subset_size, idxs + cum, W, INSERT_VALUES));
51025d06dbeSstefano_zampini         cum += subset_size;
51125d06dbeSstefano_zampini       }
5129566063dSJacob Faibussowitsch       PetscCall(MatAssemblyBegin(T, MAT_FINAL_ASSEMBLY));
5139566063dSJacob Faibussowitsch       PetscCall(MatAssemblyEnd(T, MAT_FINAL_ASSEMBLY));
5149566063dSJacob Faibussowitsch       PetscCall(MatMatTransposeMult(T, fetidpmat_ctx->B_delta, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &BD1));
5159566063dSJacob Faibussowitsch       PetscCall(MatMatMult(fetidpmat_ctx->B_delta, BD1, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &BD2));
5169566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&T));
5179566063dSJacob Faibussowitsch       PetscCall(PetscFree3(W, pivots, Bwork));
518*48a46eb9SPierre Jolivet       if (sub_schurs->is_Ej_all) PetscCall(ISRestoreIndices(sub_schurs->is_Ej_all, &idxs));
51925d06dbeSstefano_zampini     }
520674ae819SStefano Zampini   }
5219566063dSJacob Faibussowitsch   PetscCall(PetscFree(scaling_factors));
5229566063dSJacob Faibussowitsch   PetscCall(PetscFree(cols_B_delta));
523674ae819SStefano Zampini 
524457d4a33Sstefano_zampini   /* Layout of multipliers */
5259566063dSJacob Faibussowitsch   PetscCall(PetscLayoutCreate(comm, &llay));
5269566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetBlockSize(llay, 1));
5279566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetSize(llay, fetidpmat_ctx->n_lambda));
5289566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(llay));
5299566063dSJacob Faibussowitsch   PetscCall(PetscLayoutGetLocalSize(llay, &fetidpmat_ctx->n));
530457d4a33Sstefano_zampini 
531457d4a33Sstefano_zampini   /* Local work vector of multipliers */
5329566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_SELF, &fetidpmat_ctx->lambda_local));
5339566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(fetidpmat_ctx->lambda_local, n_local_lambda, n_local_lambda));
5349566063dSJacob Faibussowitsch   PetscCall(VecSetType(fetidpmat_ctx->lambda_local, VECSEQ));
535457d4a33Sstefano_zampini 
53625d06dbeSstefano_zampini   if (BD2) {
53725d06dbeSstefano_zampini     ISLocalToGlobalMapping l2g;
53825d06dbeSstefano_zampini     Mat                    T, TA, *pT;
53925d06dbeSstefano_zampini     IS                     is;
54025d06dbeSstefano_zampini     PetscInt               nl, N;
54125d06dbeSstefano_zampini     BDdelta_DN             ctx;
54225d06dbeSstefano_zampini 
5439566063dSJacob Faibussowitsch     PetscCall(PetscLayoutGetLocalSize(llay, &nl));
5449566063dSJacob Faibussowitsch     PetscCall(PetscLayoutGetSize(llay, &N));
5459566063dSJacob Faibussowitsch     PetscCall(MatCreate(comm, &T));
5469566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(T, nl, nl, N, N));
5479566063dSJacob Faibussowitsch     PetscCall(MatSetType(T, MATIS));
5489566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingCreate(comm, 1, n_local_lambda, l2g_indices, PETSC_COPY_VALUES, &l2g));
5499566063dSJacob Faibussowitsch     PetscCall(MatSetLocalToGlobalMapping(T, l2g, l2g));
5509566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingDestroy(&l2g));
5519566063dSJacob Faibussowitsch     PetscCall(MatISSetLocalMat(T, BD2));
5529566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(T, MAT_FINAL_ASSEMBLY));
5539566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(T, MAT_FINAL_ASSEMBLY));
5549566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&BD2));
5559566063dSJacob Faibussowitsch     PetscCall(MatConvert(T, MATAIJ, MAT_INITIAL_MATRIX, &TA));
5569566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&T));
5579566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(comm, n_local_lambda, l2g_indices, PETSC_USE_POINTER, &is));
5589566063dSJacob Faibussowitsch     PetscCall(MatCreateSubMatrices(TA, 1, &is, &is, MAT_INITIAL_MATRIX, &pT));
5599566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&TA));
5609566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&is));
56125d06dbeSstefano_zampini     BD2 = pT[0];
5629566063dSJacob Faibussowitsch     PetscCall(PetscFree(pT));
56325d06dbeSstefano_zampini 
56425d06dbeSstefano_zampini     /* B_Ddelta for non-redundant multipliers with deluxe scaling */
5659566063dSJacob Faibussowitsch     PetscCall(PetscNew(&ctx));
5669566063dSJacob Faibussowitsch     PetscCall(MatSetType(fetidpmat_ctx->B_Ddelta, MATSHELL));
5679566063dSJacob Faibussowitsch     PetscCall(MatShellSetContext(fetidpmat_ctx->B_Ddelta, ctx));
5689566063dSJacob Faibussowitsch     PetscCall(MatShellSetOperation(fetidpmat_ctx->B_Ddelta, MATOP_MULT, (void (*)(void))MatMult_BDdelta_deluxe_nonred));
5699566063dSJacob Faibussowitsch     PetscCall(MatShellSetOperation(fetidpmat_ctx->B_Ddelta, MATOP_MULT_TRANSPOSE, (void (*)(void))MatMultTranspose_BDdelta_deluxe_nonred));
5709566063dSJacob Faibussowitsch     PetscCall(MatShellSetOperation(fetidpmat_ctx->B_Ddelta, MATOP_DESTROY, (void (*)(void))MatDestroy_BDdelta_deluxe_nonred));
5719566063dSJacob Faibussowitsch     PetscCall(MatSetUp(fetidpmat_ctx->B_Ddelta));
57225d06dbeSstefano_zampini 
5739566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)BD1));
57425d06dbeSstefano_zampini     ctx->BD = BD1;
5759566063dSJacob Faibussowitsch     PetscCall(KSPCreate(PETSC_COMM_SELF, &ctx->kBD));
5769566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(ctx->kBD, BD2, BD2));
5779566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(fetidpmat_ctx->lambda_local, &ctx->work));
57825d06dbeSstefano_zampini     fetidpmat_ctx->deluxe_nonred = PETSC_TRUE;
57925d06dbeSstefano_zampini   }
5809566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&BD1));
5819566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&BD2));
58225d06dbeSstefano_zampini 
58325d06dbeSstefano_zampini   /* fetidpmat sizes */
58425d06dbeSstefano_zampini   fetidpmat_ctx->n += nPgl;
58525d06dbeSstefano_zampini   fetidpmat_ctx->N = fetidpmat_ctx->n_lambda + nPg;
58625d06dbeSstefano_zampini 
587457d4a33Sstefano_zampini   /* Global vector for FETI-DP linear system */
5889566063dSJacob Faibussowitsch   PetscCall(VecCreate(comm, &fetidp_global));
5899566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(fetidp_global, fetidpmat_ctx->n, fetidpmat_ctx->N));
5909566063dSJacob Faibussowitsch   PetscCall(VecSetType(fetidp_global, VECMPI));
5919566063dSJacob Faibussowitsch   PetscCall(VecSetUp(fetidp_global));
592457d4a33Sstefano_zampini 
5939eec4de8Sstefano_zampini   /* Decide layout for fetidp dofs: if it is a saddle point problem
5949eec4de8Sstefano_zampini      pressure is ordered first in the local part of the global vector
595457d4a33Sstefano_zampini      of the FETI-DP linear system */
596457d4a33Sstefano_zampini   if (nPg) {
597e1214c54Sstefano_zampini     Vec             v;
598af140850Sstefano_zampini     IS              IS_l2g_p, ais;
599457d4a33Sstefano_zampini     PetscLayout     alay;
600457d4a33Sstefano_zampini     const PetscInt *idxs, *pranges, *aranges, *lranges;
601af140850Sstefano_zampini     PetscInt       *l2g_indices_p, rst;
602e1214c54Sstefano_zampini     PetscMPIInt     rank;
603457d4a33Sstefano_zampini 
6049566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nPl, &l2g_indices_p));
6059566063dSJacob Faibussowitsch     PetscCall(VecGetLayout(fetidp_global, &alay));
6069566063dSJacob Faibussowitsch     PetscCall(PetscLayoutGetRanges(alay, &aranges));
6079566063dSJacob Faibussowitsch     PetscCall(PetscLayoutGetRanges(play, &pranges));
6089566063dSJacob Faibussowitsch     PetscCall(PetscLayoutGetRanges(llay, &lranges));
609e1214c54Sstefano_zampini 
6109566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)fetidp_global), &rank));
6119566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PetscObjectComm((PetscObject)fetidp_global), pranges[rank + 1] - pranges[rank], aranges[rank], 1, &fetidpmat_ctx->pressure));
6129566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)fetidpmat_ctx->pressure, "F_P"));
6139566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PetscObjectComm((PetscObject)fetidp_global), lranges[rank + 1] - lranges[rank], aranges[rank] + pranges[rank + 1] - pranges[rank], 1, &fetidpmat_ctx->lagrange));
6149566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)fetidpmat_ctx->lagrange, "F_L"));
6159566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingGetIndices(l2gmap_p, &idxs));
616af140850Sstefano_zampini     /* shift local to global indices for pressure */
617457d4a33Sstefano_zampini     for (i = 0; i < nPl; i++) {
618131c27b5Sprj-       PetscMPIInt owner;
619457d4a33Sstefano_zampini 
6209566063dSJacob Faibussowitsch       PetscCall(PetscLayoutFindOwner(play, idxs[i], &owner));
621457d4a33Sstefano_zampini       l2g_indices_p[i] = idxs[i] - pranges[owner] + aranges[owner];
622457d4a33Sstefano_zampini     }
6239566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingRestoreIndices(l2gmap_p, &idxs));
6249566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(comm, nPl, l2g_indices_p, PETSC_OWN_POINTER, &IS_l2g_p));
625af140850Sstefano_zampini 
626e1214c54Sstefano_zampini     /* local to global scatter for pressure */
6279566063dSJacob Faibussowitsch     PetscCall(VecScatterCreate(fetidpmat_ctx->vP, NULL, fetidp_global, IS_l2g_p, &fetidpmat_ctx->l2g_p));
6289566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&IS_l2g_p));
629457d4a33Sstefano_zampini 
630e1214c54Sstefano_zampini     /* scatter for lagrange multipliers only */
6319566063dSJacob Faibussowitsch     PetscCall(VecCreate(comm, &v));
6329566063dSJacob Faibussowitsch     PetscCall(VecSetType(v, VECSTANDARD));
6339566063dSJacob Faibussowitsch     PetscCall(VecSetLayout(v, llay));
6349566063dSJacob Faibussowitsch     PetscCall(VecSetUp(v));
6359566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(comm, n_local_lambda, l2g_indices, PETSC_COPY_VALUES, &ais));
6369566063dSJacob Faibussowitsch     PetscCall(VecScatterCreate(fetidpmat_ctx->lambda_local, NULL, v, ais, &fetidpmat_ctx->l2g_lambda_only));
6379566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&ais));
6389566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&v));
639e1214c54Sstefano_zampini 
640af140850Sstefano_zampini     /* shift local to global indices for multipliers */
641457d4a33Sstefano_zampini     for (i = 0; i < n_local_lambda; i++) {
642131c27b5Sprj-       PetscInt    ps;
643131c27b5Sprj-       PetscMPIInt owner;
644457d4a33Sstefano_zampini 
6459566063dSJacob Faibussowitsch       PetscCall(PetscLayoutFindOwner(llay, l2g_indices[i], &owner));
646457d4a33Sstefano_zampini       ps             = pranges[owner + 1] - pranges[owner];
647457d4a33Sstefano_zampini       l2g_indices[i] = l2g_indices[i] - lranges[owner] + aranges[owner] + ps;
648457d4a33Sstefano_zampini     }
649457d4a33Sstefano_zampini 
650e1214c54Sstefano_zampini     /* scatter from alldofs to pressures global fetidp vector */
6519566063dSJacob Faibussowitsch     PetscCall(PetscLayoutGetRange(alay, &rst, NULL));
6529566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(comm, nPgl, rst, 1, &ais));
6539566063dSJacob Faibussowitsch     PetscCall(VecScatterCreate(pcis->vec1_global, pP, fetidp_global, ais, &fetidpmat_ctx->g2g_p));
6549566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&ais));
655457d4a33Sstefano_zampini   }
6569566063dSJacob Faibussowitsch   PetscCall(PetscLayoutDestroy(&llay));
6579566063dSJacob Faibussowitsch   PetscCall(PetscLayoutDestroy(&play));
6589566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(comm, n_local_lambda, l2g_indices, PETSC_OWN_POINTER, &IS_l2g_lambda));
659a1c0d0daSstefano_zampini 
6609cc7774eSstefano_zampini   /* scatter from local to global multipliers */
6619566063dSJacob Faibussowitsch   PetscCall(VecScatterCreate(fetidpmat_ctx->lambda_local, NULL, fetidp_global, IS_l2g_lambda, &fetidpmat_ctx->l2g_lambda));
6629566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&IS_l2g_lambda));
6639566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingDestroy(&l2gmap_p));
6649566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&fetidp_global));
665457d4a33Sstefano_zampini 
666a1c0d0daSstefano_zampini   /* Create some work vectors needed by fetidp */
6679566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(pcis->vec1_B, &fetidpmat_ctx->temp_solution_B));
6689566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(pcis->vec1_D, &fetidpmat_ctx->temp_solution_D));
669674ae819SStefano Zampini   PetscFunctionReturn(0);
670674ae819SStefano Zampini }
671674ae819SStefano Zampini 
6729371c9d4SSatish Balay PetscErrorCode PCBDDCSetupFETIDPPCContext(Mat fetimat, FETIDPPC_ctx fetidppc_ctx) {
673674ae819SStefano Zampini   FETIDPMat_ctx mat_ctx;
67425d06dbeSstefano_zampini   PC_BDDC      *pcbddc = (PC_BDDC *)fetidppc_ctx->pc->data;
67525d06dbeSstefano_zampini   PC_IS        *pcis   = (PC_IS *)fetidppc_ctx->pc->data;
676f28b6018SStefano Zampini   PetscBool     lumped = PETSC_FALSE;
677674ae819SStefano Zampini 
678674ae819SStefano Zampini   PetscFunctionBegin;
6799566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(fetimat, &mat_ctx));
680674ae819SStefano Zampini   /* get references from objects created when setting up feti mat context */
6819566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)mat_ctx->lambda_local));
682674ae819SStefano Zampini   fetidppc_ctx->lambda_local = mat_ctx->lambda_local;
6839566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)mat_ctx->B_Ddelta));
684674ae819SStefano Zampini   fetidppc_ctx->B_Ddelta = mat_ctx->B_Ddelta;
68525d06dbeSstefano_zampini   if (mat_ctx->deluxe_nonred) {
68625d06dbeSstefano_zampini     PC            pc, mpc;
68725d06dbeSstefano_zampini     BDdelta_DN    ctx;
6883ca39a21SBarry Smith     MatSolverType solver;
68925d06dbeSstefano_zampini     const char   *prefix;
69025d06dbeSstefano_zampini 
6919566063dSJacob Faibussowitsch     PetscCall(MatShellGetContext(mat_ctx->B_Ddelta, &ctx));
6929566063dSJacob Faibussowitsch     PetscCall(KSPSetType(ctx->kBD, KSPPREONLY));
6939566063dSJacob Faibussowitsch     PetscCall(KSPGetPC(ctx->kBD, &mpc));
6949566063dSJacob Faibussowitsch     PetscCall(KSPGetPC(pcbddc->ksp_D, &pc));
6959566063dSJacob Faibussowitsch     PetscCall(PCSetType(mpc, PCLU));
6969566063dSJacob Faibussowitsch     PetscCall(PCFactorGetMatSolverType(pc, (MatSolverType *)&solver));
6971baa6e33SBarry Smith     if (solver) PetscCall(PCFactorSetMatSolverType(mpc, solver));
6989566063dSJacob Faibussowitsch     PetscCall(MatGetOptionsPrefix(fetimat, &prefix));
6999566063dSJacob Faibussowitsch     PetscCall(KSPSetOptionsPrefix(ctx->kBD, prefix));
7009566063dSJacob Faibussowitsch     PetscCall(KSPAppendOptionsPrefix(ctx->kBD, "bddelta_"));
7019566063dSJacob Faibussowitsch     PetscCall(KSPSetFromOptions(ctx->kBD));
70225d06dbeSstefano_zampini   }
70325d06dbeSstefano_zampini 
704e1214c54Sstefano_zampini   if (mat_ctx->l2g_lambda_only) {
7059566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)mat_ctx->l2g_lambda_only));
706e1214c54Sstefano_zampini     fetidppc_ctx->l2g_lambda = mat_ctx->l2g_lambda_only;
707e1214c54Sstefano_zampini   } else {
7089566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)mat_ctx->l2g_lambda));
709674ae819SStefano Zampini     fetidppc_ctx->l2g_lambda = mat_ctx->l2g_lambda;
710e1214c54Sstefano_zampini   }
711f28b6018SStefano Zampini   /* Dirichlet preconditioner */
7129566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, ((PetscObject)fetimat)->prefix, "-pc_lumped", &lumped, NULL));
713f28b6018SStefano Zampini   if (!lumped) {
7149feb908dSstefano_zampini     IS        iV;
7159c2d02cdSstefano_zampini     PetscBool discrete_harmonic = PETSC_FALSE;
7169c2d02cdSstefano_zampini 
7179566063dSJacob Faibussowitsch     PetscCall(PetscObjectQuery((PetscObject)fetidppc_ctx->pc, "__KSPFETIDP_iV", (PetscObject *)&iV));
718*48a46eb9SPierre Jolivet     if (iV) PetscCall(PetscOptionsGetBool(NULL, ((PetscObject)fetimat)->prefix, "-pc_discrete_harmonic", &discrete_harmonic, NULL));
7199c2d02cdSstefano_zampini     if (discrete_harmonic) {
7209c2d02cdSstefano_zampini       KSP             sksp;
7219c2d02cdSstefano_zampini       PC              pc;
722111315fdSstefano_zampini       PCBDDCSubSchurs sub_schurs = pcbddc->sub_schurs;
7239c2d02cdSstefano_zampini       Mat             A_II, A_IB, A_BI;
724111315fdSstefano_zampini       IS              iP = NULL;
725111315fdSstefano_zampini       PetscBool       isshell, reuse = PETSC_FALSE;
7269c2d02cdSstefano_zampini       KSPType         ksptype;
7279c2d02cdSstefano_zampini       const char     *prefix;
7289c2d02cdSstefano_zampini 
7299c2d02cdSstefano_zampini       /*
7309c2d02cdSstefano_zampini         We constructs a Schur complement for
7319c2d02cdSstefano_zampini 
7329c2d02cdSstefano_zampini         | A_II A_ID |
7339c2d02cdSstefano_zampini         | A_DI A_DD |
7349c2d02cdSstefano_zampini 
7359c2d02cdSstefano_zampini         instead of
7369c2d02cdSstefano_zampini 
7379c2d02cdSstefano_zampini         | A_II  B^t_II A_ID |
7389c2d02cdSstefano_zampini         | B_II -C_II   B_ID |
7399c2d02cdSstefano_zampini         | A_DI  B^t_ID A_DD |
7409c2d02cdSstefano_zampini 
7419c2d02cdSstefano_zampini       */
742111315fdSstefano_zampini       if (sub_schurs && sub_schurs->reuse_solver) {
7439566063dSJacob Faibussowitsch         PetscCall(PetscObjectQuery((PetscObject)sub_schurs->A, "__KSPFETIDP_iP", (PetscObject *)&iP));
744111315fdSstefano_zampini         if (iP) reuse = PETSC_TRUE;
745111315fdSstefano_zampini       }
746111315fdSstefano_zampini       if (!reuse) {
747111315fdSstefano_zampini         IS       aB;
748111315fdSstefano_zampini         PetscInt nb;
7499566063dSJacob Faibussowitsch         PetscCall(ISGetLocalSize(pcis->is_B_local, &nb));
7509566063dSJacob Faibussowitsch         PetscCall(ISCreateStride(PetscObjectComm((PetscObject)pcis->A_II), nb, 0, 1, &aB));
7519566063dSJacob Faibussowitsch         PetscCall(MatCreateSubMatrix(pcis->A_II, iV, iV, MAT_INITIAL_MATRIX, &A_II));
7529566063dSJacob Faibussowitsch         PetscCall(MatCreateSubMatrix(pcis->A_IB, iV, aB, MAT_INITIAL_MATRIX, &A_IB));
7539566063dSJacob Faibussowitsch         PetscCall(MatCreateSubMatrix(pcis->A_BI, aB, iV, MAT_INITIAL_MATRIX, &A_BI));
7549566063dSJacob Faibussowitsch         PetscCall(ISDestroy(&aB));
755111315fdSstefano_zampini       } else {
7569566063dSJacob Faibussowitsch         PetscCall(MatCreateSubMatrix(sub_schurs->A, pcis->is_I_local, pcis->is_B_local, MAT_INITIAL_MATRIX, &A_IB));
7579566063dSJacob Faibussowitsch         PetscCall(MatCreateSubMatrix(sub_schurs->A, pcis->is_B_local, pcis->is_I_local, MAT_INITIAL_MATRIX, &A_BI));
7589566063dSJacob Faibussowitsch         PetscCall(PetscObjectReference((PetscObject)pcis->A_II));
759111315fdSstefano_zampini         A_II = pcis->A_II;
760111315fdSstefano_zampini       }
7619566063dSJacob Faibussowitsch       PetscCall(MatCreateSchurComplement(A_II, A_II, A_IB, A_BI, pcis->A_BB, &fetidppc_ctx->S_j));
7629c2d02cdSstefano_zampini 
7639c2d02cdSstefano_zampini       /* propagate settings of solver */
7649566063dSJacob Faibussowitsch       PetscCall(MatSchurComplementGetKSP(fetidppc_ctx->S_j, &sksp));
7659566063dSJacob Faibussowitsch       PetscCall(KSPGetType(pcis->ksp_D, &ksptype));
7669566063dSJacob Faibussowitsch       PetscCall(KSPSetType(sksp, ksptype));
7679566063dSJacob Faibussowitsch       PetscCall(KSPGetPC(pcis->ksp_D, &pc));
7689566063dSJacob Faibussowitsch       PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCSHELL, &isshell));
7695334bea6Sstefano_zampini       if (!isshell) {
7703ca39a21SBarry Smith         MatSolverType solver;
7715334bea6Sstefano_zampini         PCType        pctype;
7725334bea6Sstefano_zampini 
7739566063dSJacob Faibussowitsch         PetscCall(PCGetType(pc, &pctype));
7749566063dSJacob Faibussowitsch         PetscCall(PCFactorGetMatSolverType(pc, (MatSolverType *)&solver));
7759566063dSJacob Faibussowitsch         PetscCall(KSPGetPC(sksp, &pc));
7769566063dSJacob Faibussowitsch         PetscCall(PCSetType(pc, pctype));
7771baa6e33SBarry Smith         if (solver) PetscCall(PCFactorSetMatSolverType(pc, solver));
7785334bea6Sstefano_zampini       } else {
7799566063dSJacob Faibussowitsch         PetscCall(KSPGetPC(sksp, &pc));
7809566063dSJacob Faibussowitsch         PetscCall(PCSetType(pc, PCLU));
7815334bea6Sstefano_zampini       }
7829566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&A_II));
7839566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&A_IB));
7849566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&A_BI));
7859566063dSJacob Faibussowitsch       PetscCall(MatGetOptionsPrefix(fetimat, &prefix));
7869566063dSJacob Faibussowitsch       PetscCall(KSPSetOptionsPrefix(sksp, prefix));
7879566063dSJacob Faibussowitsch       PetscCall(KSPAppendOptionsPrefix(sksp, "harmonic_"));
7889566063dSJacob Faibussowitsch       PetscCall(KSPSetFromOptions(sksp));
789111315fdSstefano_zampini       if (reuse) {
7909566063dSJacob Faibussowitsch         PetscCall(KSPSetPC(sksp, sub_schurs->reuse_solver->interior_solver));
7919566063dSJacob Faibussowitsch         PetscCall(PetscObjectIncrementTabLevel((PetscObject)sub_schurs->reuse_solver->interior_solver, (PetscObject)sksp, 0));
792111315fdSstefano_zampini       }
7939c2d02cdSstefano_zampini     } else { /* default Dirichlet preconditioner is pde-harmonic */
7949566063dSJacob Faibussowitsch       PetscCall(MatCreateSchurComplement(pcis->A_II, pcis->A_II, pcis->A_IB, pcis->A_BI, pcis->A_BB, &fetidppc_ctx->S_j));
7959566063dSJacob Faibussowitsch       PetscCall(MatSchurComplementSetKSP(fetidppc_ctx->S_j, pcis->ksp_D));
7969c2d02cdSstefano_zampini     }
797f28b6018SStefano Zampini   } else {
7989566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)pcis->A_BB));
799f28b6018SStefano Zampini     fetidppc_ctx->S_j = pcis->A_BB;
800f28b6018SStefano Zampini   }
801af140850Sstefano_zampini   /* saddle-point */
802af140850Sstefano_zampini   if (mat_ctx->xPg) {
8039566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)mat_ctx->xPg));
804af140850Sstefano_zampini     fetidppc_ctx->xPg = mat_ctx->xPg;
8059566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)mat_ctx->yPg));
806af140850Sstefano_zampini     fetidppc_ctx->yPg = mat_ctx->yPg;
807af140850Sstefano_zampini   }
808674ae819SStefano Zampini   PetscFunctionReturn(0);
809674ae819SStefano Zampini }
810674ae819SStefano Zampini 
8119371c9d4SSatish Balay PetscErrorCode FETIDPMatMult_Kernel(Mat fetimat, Vec x, Vec y, PetscBool trans) {
812674ae819SStefano Zampini   FETIDPMat_ctx mat_ctx;
813617d11aeSStefano Zampini   PC_BDDC      *pcbddc;
814674ae819SStefano Zampini   PC_IS        *pcis;
815674ae819SStefano Zampini 
816674ae819SStefano Zampini   PetscFunctionBegin;
8179566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(fetimat, &mat_ctx));
818674ae819SStefano Zampini   pcis   = (PC_IS *)mat_ctx->pc->data;
819617d11aeSStefano Zampini   pcbddc = (PC_BDDC *)mat_ctx->pc->data;
820674ae819SStefano Zampini   /* Application of B_delta^T */
8219566063dSJacob Faibussowitsch   PetscCall(VecSet(pcis->vec1_B, 0.));
8229566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(mat_ctx->l2g_lambda, x, mat_ctx->lambda_local, INSERT_VALUES, SCATTER_REVERSE));
8239566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(mat_ctx->l2g_lambda, x, mat_ctx->lambda_local, INSERT_VALUES, SCATTER_REVERSE));
8249566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(mat_ctx->B_delta, mat_ctx->lambda_local, pcis->vec1_B));
825af140850Sstefano_zampini 
826af140850Sstefano_zampini   /* Add contribution from saddle point */
827af140850Sstefano_zampini   if (mat_ctx->l2g_p) {
8289566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(mat_ctx->l2g_p, x, mat_ctx->vP, INSERT_VALUES, SCATTER_REVERSE));
8299566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(mat_ctx->l2g_p, x, mat_ctx->vP, INSERT_VALUES, SCATTER_REVERSE));
830af140850Sstefano_zampini     if (pcbddc->switch_static) {
831ef425aeaSstefano_zampini       if (trans) {
8329566063dSJacob Faibussowitsch         PetscCall(MatMultTranspose(mat_ctx->B_BI, mat_ctx->vP, pcis->vec1_D));
833ef425aeaSstefano_zampini       } else {
8349566063dSJacob Faibussowitsch         PetscCall(MatMult(mat_ctx->Bt_BI, mat_ctx->vP, pcis->vec1_D));
835af140850Sstefano_zampini       }
836ef425aeaSstefano_zampini     }
837ef425aeaSstefano_zampini     if (trans) {
8389566063dSJacob Faibussowitsch       PetscCall(MatMultTransposeAdd(mat_ctx->B_BB, mat_ctx->vP, pcis->vec1_B, pcis->vec1_B));
839ef425aeaSstefano_zampini     } else {
8409566063dSJacob Faibussowitsch       PetscCall(MatMultAdd(mat_ctx->Bt_BB, mat_ctx->vP, pcis->vec1_B, pcis->vec1_B));
841ef425aeaSstefano_zampini     }
842af140850Sstefano_zampini   } else {
8431baa6e33SBarry Smith     if (pcbddc->switch_static) PetscCall(VecSet(pcis->vec1_D, 0.0));
844af140850Sstefano_zampini   }
845af140850Sstefano_zampini   /* Application of \widetilde{S}^-1 */
8469566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(pcbddc->benign_p0, pcbddc->benign_n));
8479566063dSJacob Faibussowitsch   PetscCall(PCBDDCApplyInterfacePreconditioner(mat_ctx->pc, trans));
8489566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(pcbddc->benign_p0, pcbddc->benign_n));
8499566063dSJacob Faibussowitsch   PetscCall(VecSet(y, 0.0));
850674ae819SStefano Zampini   /* Application of B_delta */
8519566063dSJacob Faibussowitsch   PetscCall(MatMult(mat_ctx->B_delta, pcis->vec1_B, mat_ctx->lambda_local));
852af140850Sstefano_zampini   /* Contribution from boundary pressures */
853af140850Sstefano_zampini   if (mat_ctx->C) {
854af140850Sstefano_zampini     const PetscScalar *lx;
855af140850Sstefano_zampini     PetscScalar       *ly;
856af140850Sstefano_zampini 
85715579a77SStefano Zampini     /* pressure ordered first in the local part of x and y */
8589566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(x, &lx));
8599566063dSJacob Faibussowitsch     PetscCall(VecGetArray(y, &ly));
8609566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(mat_ctx->xPg, lx));
8619566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(mat_ctx->yPg, ly));
862ef425aeaSstefano_zampini     if (trans) {
8639566063dSJacob Faibussowitsch       PetscCall(MatMultTranspose(mat_ctx->C, mat_ctx->xPg, mat_ctx->yPg));
864ef425aeaSstefano_zampini     } else {
8659566063dSJacob Faibussowitsch       PetscCall(MatMult(mat_ctx->C, mat_ctx->xPg, mat_ctx->yPg));
866ef425aeaSstefano_zampini     }
8679566063dSJacob Faibussowitsch     PetscCall(VecResetArray(mat_ctx->xPg));
8689566063dSJacob Faibussowitsch     PetscCall(VecResetArray(mat_ctx->yPg));
8699566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(x, &lx));
8709566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(y, &ly));
871af140850Sstefano_zampini   }
872af140850Sstefano_zampini   /* Add contribution from saddle point */
873af140850Sstefano_zampini   if (mat_ctx->l2g_p) {
874ef425aeaSstefano_zampini     if (trans) {
8759566063dSJacob Faibussowitsch       PetscCall(MatMultTranspose(mat_ctx->Bt_BB, pcis->vec1_B, mat_ctx->vP));
876ef425aeaSstefano_zampini     } else {
8779566063dSJacob Faibussowitsch       PetscCall(MatMult(mat_ctx->B_BB, pcis->vec1_B, mat_ctx->vP));
878ef425aeaSstefano_zampini     }
879af140850Sstefano_zampini     if (pcbddc->switch_static) {
880ef425aeaSstefano_zampini       if (trans) {
8819566063dSJacob Faibussowitsch         PetscCall(MatMultTransposeAdd(mat_ctx->Bt_BI, pcis->vec1_D, mat_ctx->vP, mat_ctx->vP));
882ef425aeaSstefano_zampini       } else {
8839566063dSJacob Faibussowitsch         PetscCall(MatMultAdd(mat_ctx->B_BI, pcis->vec1_D, mat_ctx->vP, mat_ctx->vP));
884af140850Sstefano_zampini       }
885ef425aeaSstefano_zampini     }
8869566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(mat_ctx->l2g_p, mat_ctx->vP, y, ADD_VALUES, SCATTER_FORWARD));
8879566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(mat_ctx->l2g_p, mat_ctx->vP, y, ADD_VALUES, SCATTER_FORWARD));
888af140850Sstefano_zampini   }
8899566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(mat_ctx->l2g_lambda, mat_ctx->lambda_local, y, ADD_VALUES, SCATTER_FORWARD));
8909566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(mat_ctx->l2g_lambda, mat_ctx->lambda_local, y, ADD_VALUES, SCATTER_FORWARD));
891674ae819SStefano Zampini   PetscFunctionReturn(0);
892674ae819SStefano Zampini }
893674ae819SStefano Zampini 
8949371c9d4SSatish Balay PetscErrorCode FETIDPMatMult(Mat fetimat, Vec x, Vec y) {
895edf7251bSStefano Zampini   PetscFunctionBegin;
8969566063dSJacob Faibussowitsch   PetscCall(FETIDPMatMult_Kernel(fetimat, x, y, PETSC_FALSE));
897ef425aeaSstefano_zampini   PetscFunctionReturn(0);
898ef425aeaSstefano_zampini }
899ef425aeaSstefano_zampini 
9009371c9d4SSatish Balay PetscErrorCode FETIDPMatMultTranspose(Mat fetimat, Vec x, Vec y) {
901ef425aeaSstefano_zampini   PetscFunctionBegin;
9029566063dSJacob Faibussowitsch   PetscCall(FETIDPMatMult_Kernel(fetimat, x, y, PETSC_TRUE));
903edf7251bSStefano Zampini   PetscFunctionReturn(0);
904edf7251bSStefano Zampini }
905edf7251bSStefano Zampini 
9069371c9d4SSatish Balay PetscErrorCode FETIDPPCApply_Kernel(PC fetipc, Vec x, Vec y, PetscBool trans) {
907674ae819SStefano Zampini   FETIDPPC_ctx pc_ctx;
908674ae819SStefano Zampini   PC_IS       *pcis;
909674ae819SStefano Zampini 
910674ae819SStefano Zampini   PetscFunctionBegin;
9119566063dSJacob Faibussowitsch   PetscCall(PCShellGetContext(fetipc, &pc_ctx));
912674ae819SStefano Zampini   pcis = (PC_IS *)pc_ctx->pc->data;
913674ae819SStefano Zampini   /* Application of B_Ddelta^T */
9149566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(pc_ctx->l2g_lambda, x, pc_ctx->lambda_local, INSERT_VALUES, SCATTER_REVERSE));
9159566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(pc_ctx->l2g_lambda, x, pc_ctx->lambda_local, INSERT_VALUES, SCATTER_REVERSE));
9169566063dSJacob Faibussowitsch   PetscCall(VecSet(pcis->vec2_B, 0.0));
9179566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(pc_ctx->B_Ddelta, pc_ctx->lambda_local, pcis->vec2_B));
918ed6c3d69SStefano Zampini   /* Application of local Schur complement */
91916597391Sstefano_zampini   if (trans) {
9209566063dSJacob Faibussowitsch     PetscCall(MatMultTranspose(pc_ctx->S_j, pcis->vec2_B, pcis->vec1_B));
92116597391Sstefano_zampini   } else {
9229566063dSJacob Faibussowitsch     PetscCall(MatMult(pc_ctx->S_j, pcis->vec2_B, pcis->vec1_B));
92316597391Sstefano_zampini   }
924edf7251bSStefano Zampini   /* Application of B_Ddelta */
9259566063dSJacob Faibussowitsch   PetscCall(MatMult(pc_ctx->B_Ddelta, pcis->vec1_B, pc_ctx->lambda_local));
9269566063dSJacob Faibussowitsch   PetscCall(VecSet(y, 0.0));
9279566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(pc_ctx->l2g_lambda, pc_ctx->lambda_local, y, ADD_VALUES, SCATTER_FORWARD));
9289566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(pc_ctx->l2g_lambda, pc_ctx->lambda_local, y, ADD_VALUES, SCATTER_FORWARD));
929edf7251bSStefano Zampini   PetscFunctionReturn(0);
930edf7251bSStefano Zampini }
931edf7251bSStefano Zampini 
9329371c9d4SSatish Balay PetscErrorCode FETIDPPCApply(PC pc, Vec x, Vec y) {
933edf7251bSStefano Zampini   PetscFunctionBegin;
9349566063dSJacob Faibussowitsch   PetscCall(FETIDPPCApply_Kernel(pc, x, y, PETSC_FALSE));
93516597391Sstefano_zampini   PetscFunctionReturn(0);
93616597391Sstefano_zampini }
93716597391Sstefano_zampini 
9389371c9d4SSatish Balay PetscErrorCode FETIDPPCApplyTranspose(PC pc, Vec x, Vec y) {
93916597391Sstefano_zampini   PetscFunctionBegin;
9409566063dSJacob Faibussowitsch   PetscCall(FETIDPPCApply_Kernel(pc, x, y, PETSC_TRUE));
941674ae819SStefano Zampini   PetscFunctionReturn(0);
942674ae819SStefano Zampini }
943c45b8d2dSstefano_zampini 
9449371c9d4SSatish Balay PetscErrorCode FETIDPPCView(PC pc, PetscViewer viewer) {
945c45b8d2dSstefano_zampini   FETIDPPC_ctx pc_ctx;
946c45b8d2dSstefano_zampini   PetscBool    iascii;
947c45b8d2dSstefano_zampini   PetscViewer  sviewer;
948c45b8d2dSstefano_zampini 
949c45b8d2dSstefano_zampini   PetscFunctionBegin;
9509566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
951c45b8d2dSstefano_zampini   if (iascii) {
952c45b8d2dSstefano_zampini     PetscMPIInt rank;
95325d06dbeSstefano_zampini     PetscBool   isschur, isshell;
954c45b8d2dSstefano_zampini 
9559566063dSJacob Faibussowitsch     PetscCall(PCShellGetContext(pc, &pc_ctx));
9569566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank));
9579566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)pc_ctx->S_j, MATSCHURCOMPLEMENT, &isschur));
958c45b8d2dSstefano_zampini     if (isschur) {
9599566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "  Dirichlet preconditioner (just from rank 0)\n"));
960c45b8d2dSstefano_zampini     } else {
9619566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "  Lumped preconditioner (just from rank 0)\n"));
962c45b8d2dSstefano_zampini     }
9639566063dSJacob Faibussowitsch     PetscCall(PetscViewerGetSubViewer(viewer, PetscObjectComm((PetscObject)pc_ctx->S_j), &sviewer));
964dd400576SPatrick Sanan     if (rank == 0) {
9659566063dSJacob Faibussowitsch       PetscCall(PetscViewerPushFormat(sviewer, PETSC_VIEWER_ASCII_INFO));
9669566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPushTab(sviewer));
9679566063dSJacob Faibussowitsch       PetscCall(MatView(pc_ctx->S_j, sviewer));
9689566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPopTab(sviewer));
9699566063dSJacob Faibussowitsch       PetscCall(PetscViewerPopFormat(sviewer));
970c45b8d2dSstefano_zampini     }
9719566063dSJacob Faibussowitsch     PetscCall(PetscViewerRestoreSubViewer(viewer, PetscObjectComm((PetscObject)pc_ctx->S_j), &sviewer));
9729566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)pc_ctx->B_Ddelta, MATSHELL, &isshell));
97325d06dbeSstefano_zampini     if (isshell) {
97425d06dbeSstefano_zampini       BDdelta_DN ctx;
9759566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "  FETI-DP BDdelta: DB^t * (B D^-1 B^t)^-1 for deluxe scaling (just from rank 0)\n"));
9769566063dSJacob Faibussowitsch       PetscCall(MatShellGetContext(pc_ctx->B_Ddelta, &ctx));
9779566063dSJacob Faibussowitsch       PetscCall(PetscViewerGetSubViewer(viewer, PetscObjectComm((PetscObject)pc_ctx->S_j), &sviewer));
978dd400576SPatrick Sanan       if (rank == 0) {
97915579a77SStefano Zampini         PetscInt tl;
98015579a77SStefano Zampini 
9819566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIGetTab(sviewer, &tl));
9829566063dSJacob Faibussowitsch         PetscCall(PetscObjectSetTabLevel((PetscObject)ctx->kBD, tl));
9839566063dSJacob Faibussowitsch         PetscCall(KSPView(ctx->kBD, sviewer));
9849566063dSJacob Faibussowitsch         PetscCall(PetscViewerPushFormat(sviewer, PETSC_VIEWER_ASCII_INFO));
9859566063dSJacob Faibussowitsch         PetscCall(MatView(ctx->BD, sviewer));
9869566063dSJacob Faibussowitsch         PetscCall(PetscViewerPopFormat(sviewer));
98725d06dbeSstefano_zampini       }
9889566063dSJacob Faibussowitsch       PetscCall(PetscViewerRestoreSubViewer(viewer, PetscObjectComm((PetscObject)pc_ctx->S_j), &sviewer));
989e5afcf28SBarry Smith     }
9909566063dSJacob Faibussowitsch     PetscCall(PetscViewerFlush(viewer));
991c45b8d2dSstefano_zampini   }
992c45b8d2dSstefano_zampini   PetscFunctionReturn(0);
993c45b8d2dSstefano_zampini }
994