xref: /petsc/src/ksp/pc/impls/bddc/bddcfetidp.c (revision 792fecdfe9134cce4d631112660ddd34f063bc17)
15e5bbd0aSStefano Zampini #include <petsc/private/pcbddcimpl.h>
25e5bbd0aSStefano Zampini #include <petsc/private/pcbddcprivateimpl.h>
325d06dbeSstefano_zampini #include <petscblaslapack.h>
425d06dbeSstefano_zampini 
525d06dbeSstefano_zampini static PetscErrorCode MatMult_BDdelta_deluxe_nonred(Mat A, Vec x, Vec y)
625d06dbeSstefano_zampini {
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() */
1425d06dbeSstefano_zampini   PetscFunctionReturn(0);
1525d06dbeSstefano_zampini }
1625d06dbeSstefano_zampini 
1725d06dbeSstefano_zampini static PetscErrorCode MatMultTranspose_BDdelta_deluxe_nonred(Mat A, Vec x, Vec y)
1825d06dbeSstefano_zampini {
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));
2625d06dbeSstefano_zampini   PetscFunctionReturn(0);
2725d06dbeSstefano_zampini }
2825d06dbeSstefano_zampini 
2925d06dbeSstefano_zampini static PetscErrorCode MatDestroy_BDdelta_deluxe_nonred(Mat A)
3025d06dbeSstefano_zampini {
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));
3925d06dbeSstefano_zampini   PetscFunctionReturn(0);
4025d06dbeSstefano_zampini }
4125d06dbeSstefano_zampini 
42674ae819SStefano Zampini PetscErrorCode PCBDDCCreateFETIDPMatContext(PC pc, FETIDPMat_ctx *fetidpmat_ctx)
43674ae819SStefano Zampini {
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;
52674ae819SStefano Zampini   PetscFunctionReturn(0);
53674ae819SStefano Zampini }
54674ae819SStefano Zampini 
55674ae819SStefano Zampini PetscErrorCode PCBDDCCreateFETIDPPCContext(PC pc, FETIDPPC_ctx *fetidppc_ctx)
56674ae819SStefano Zampini {
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;
65674ae819SStefano Zampini   PetscFunctionReturn(0);
66674ae819SStefano Zampini }
67674ae819SStefano Zampini 
68674ae819SStefano Zampini PetscErrorCode PCBDDCDestroyFETIDPMat(Mat A)
69674ae819SStefano Zampini {
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));
96674ae819SStefano Zampini   PetscFunctionReturn(0);
97674ae819SStefano Zampini }
98674ae819SStefano Zampini 
99674ae819SStefano Zampini PetscErrorCode PCBDDCDestroyFETIDPPC(PC pc)
100674ae819SStefano Zampini {
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));
113674ae819SStefano Zampini   PetscFunctionReturn(0);
114674ae819SStefano Zampini }
115674ae819SStefano Zampini 
116674ae819SStefano Zampini PetscErrorCode PCBDDCSetupFETIDPMatContext(FETIDPMat_ctx fetidpmat_ctx)
117674ae819SStefano Zampini {
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;
139457d4a33Sstefano_zampini   PetscLayout    llay;
140a1c0d0daSstefano_zampini 
141457d4a33Sstefano_zampini   /* saddlepoint */
142457d4a33Sstefano_zampini   ISLocalToGlobalMapping l2gmap_p;
143457d4a33Sstefano_zampini   PetscLayout            play;
144457d4a33Sstefano_zampini   IS                     gP,pP;
145457d4a33Sstefano_zampini   PetscInt               nPl,nPg,nPgl;
146674ae819SStefano Zampini 
147674ae819SStefano Zampini   PetscFunctionBegin;
1489566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)(fetidpmat_ctx->pc),&comm));
1499566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm,&rank));
1509566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm,&size));
151674ae819SStefano Zampini 
152457d4a33Sstefano_zampini   /* saddlepoint */
153457d4a33Sstefano_zampini   nPl      = 0;
154457d4a33Sstefano_zampini   nPg      = 0;
155457d4a33Sstefano_zampini   nPgl     = 0;
156457d4a33Sstefano_zampini   gP       = NULL;
157457d4a33Sstefano_zampini   pP       = NULL;
158457d4a33Sstefano_zampini   l2gmap_p = NULL;
159457d4a33Sstefano_zampini   play     = NULL;
1609566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_pP",(PetscObject*)&pP));
161022d8d2bSstefano_zampini   if (pP) { /* saddle point */
162457d4a33Sstefano_zampini     /* subdomain pressures in global numbering */
1639566063dSJacob Faibussowitsch     PetscCall(PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_gP",(PetscObject*)&gP));
16428b400f6SJacob Faibussowitsch     PetscCheck(gP,PETSC_COMM_SELF,PETSC_ERR_PLIB,"gP not present");
1659566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(gP,&nPl));
1669566063dSJacob Faibussowitsch     PetscCall(VecCreate(PETSC_COMM_SELF,&fetidpmat_ctx->vP));
1679566063dSJacob Faibussowitsch     PetscCall(VecSetSizes(fetidpmat_ctx->vP,nPl,nPl));
1689566063dSJacob Faibussowitsch     PetscCall(VecSetType(fetidpmat_ctx->vP,VECSTANDARD));
1699566063dSJacob Faibussowitsch     PetscCall(VecSetUp(fetidpmat_ctx->vP));
170457d4a33Sstefano_zampini 
171e1214c54Sstefano_zampini     /* pressure matrix */
1729566063dSJacob Faibussowitsch     PetscCall(PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_C",(PetscObject*)&fetidpmat_ctx->C));
173e1214c54Sstefano_zampini     if (!fetidpmat_ctx->C) { /* null pressure block, compute layout and global numbering for pressures */
174457d4a33Sstefano_zampini       IS Pg;
175457d4a33Sstefano_zampini 
1769566063dSJacob Faibussowitsch       PetscCall(ISRenumber(gP,NULL,&nPg,&Pg));
1779566063dSJacob Faibussowitsch       PetscCall(ISLocalToGlobalMappingCreateIS(Pg,&l2gmap_p));
1789566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&Pg));
1799566063dSJacob Faibussowitsch       PetscCall(PetscLayoutCreate(comm,&play));
1809566063dSJacob Faibussowitsch       PetscCall(PetscLayoutSetBlockSize(play,1));
1819566063dSJacob Faibussowitsch       PetscCall(PetscLayoutSetSize(play,nPg));
1829566063dSJacob Faibussowitsch       PetscCall(ISGetLocalSize(pP,&nPgl));
1839566063dSJacob Faibussowitsch       PetscCall(PetscLayoutSetLocalSize(play,nPgl));
1849566063dSJacob Faibussowitsch       PetscCall(PetscLayoutSetUp(play));
185457d4a33Sstefano_zampini     } else {
1869566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)fetidpmat_ctx->C));
1879566063dSJacob Faibussowitsch       PetscCall(MatISGetLocalToGlobalMapping(fetidpmat_ctx->C,&l2gmap_p,NULL));
1889566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)l2gmap_p));
1899566063dSJacob Faibussowitsch       PetscCall(MatGetSize(fetidpmat_ctx->C,&nPg,NULL));
1909566063dSJacob Faibussowitsch       PetscCall(MatGetLocalSize(fetidpmat_ctx->C,NULL,&nPgl));
1919566063dSJacob Faibussowitsch       PetscCall(MatGetLayouts(fetidpmat_ctx->C,NULL,&llay));
1929566063dSJacob Faibussowitsch       PetscCall(PetscLayoutReference(llay,&play));
193457d4a33Sstefano_zampini     }
1949566063dSJacob Faibussowitsch     PetscCall(VecCreateMPIWithArray(comm,1,nPgl,nPg,NULL,&fetidpmat_ctx->xPg));
1959566063dSJacob Faibussowitsch     PetscCall(VecCreateMPIWithArray(comm,1,nPgl,nPg,NULL,&fetidpmat_ctx->yPg));
196457d4a33Sstefano_zampini 
197e1214c54Sstefano_zampini     /* import matrices for pressures coupling */
1989566063dSJacob Faibussowitsch     PetscCall(PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_B_BI",(PetscObject*)&fetidpmat_ctx->B_BI));
19928b400f6SJacob Faibussowitsch     PetscCheck(fetidpmat_ctx->B_BI,PETSC_COMM_SELF,PETSC_ERR_PLIB,"B_BI not present");
2009566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)fetidpmat_ctx->B_BI));
201457d4a33Sstefano_zampini 
2029566063dSJacob Faibussowitsch     PetscCall(PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_B_BB",(PetscObject*)&fetidpmat_ctx->B_BB));
20328b400f6SJacob Faibussowitsch     PetscCheck(fetidpmat_ctx->B_BB,PETSC_COMM_SELF,PETSC_ERR_PLIB,"B_BB not present");
2049566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)fetidpmat_ctx->B_BB));
205457d4a33Sstefano_zampini 
2069566063dSJacob Faibussowitsch     PetscCall(PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_Bt_BI",(PetscObject*)&fetidpmat_ctx->Bt_BI));
20728b400f6SJacob Faibussowitsch     PetscCheck(fetidpmat_ctx->Bt_BI,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Bt_BI not present");
2089566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)fetidpmat_ctx->Bt_BI));
209457d4a33Sstefano_zampini 
2109566063dSJacob Faibussowitsch     PetscCall(PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_Bt_BB",(PetscObject*)&fetidpmat_ctx->Bt_BB));
21128b400f6SJacob Faibussowitsch     PetscCheck(fetidpmat_ctx->Bt_BB,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Bt_BB not present");
2129566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)fetidpmat_ctx->Bt_BB));
2136cc1294bSstefano_zampini 
2149566063dSJacob Faibussowitsch     PetscCall(PetscObjectQuery((PetscObject)fetidpmat_ctx->pc,"__KSPFETIDP_flip" ,(PetscObject*)&fetidpmat_ctx->rhs_flip));
2151baa6e33SBarry Smith     if (fetidpmat_ctx->rhs_flip) PetscCall(PetscObjectReference((PetscObject)fetidpmat_ctx->rhs_flip));
216457d4a33Sstefano_zampini   }
217457d4a33Sstefano_zampini 
218674ae819SStefano Zampini   /* Default type of lagrange multipliers is non-redundant */
219329cd39dSStefano Zampini   fully_redundant = fetidpmat_ctx->fully_redundant;
220674ae819SStefano Zampini 
221674ae819SStefano Zampini   /* Evaluate local and global number of lagrange multipliers */
2229566063dSJacob Faibussowitsch   PetscCall(VecSet(pcis->vec1_N,0.0));
223674ae819SStefano Zampini   n_local_lambda = 0;
224674ae819SStefano Zampini   partial_sum = 0;
225674ae819SStefano Zampini   n_boundary_dofs = 0;
226674ae819SStefano Zampini   s = 0;
227a1c0d0daSstefano_zampini 
228674ae819SStefano Zampini   /* Get Vertices used to define the BDDC */
2299566063dSJacob Faibussowitsch   PetscCall(PCBDDCGraphGetCandidatesIS(pcbddc->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 
2389566063dSJacob Faibussowitsch   PetscCall(VecGetArray(pcis->vec1_N,&array));
239674ae819SStefano Zampini   for (i=0;i<pcis->n;i++) {
240674ae819SStefano Zampini     j = mat_graph->count[i]; /* RECALL: mat_graph->count[i] does not count myself */
2412d456bbaSstefano_zampini     if (j > 0) n_boundary_dofs++;
242674ae819SStefano Zampini     skip_node = PETSC_FALSE;
243674ae819SStefano Zampini     if (s < n_vertices && vertex_indices[s] == i) { /* it works for a sorted set of vertices */
244674ae819SStefano Zampini       skip_node = PETSC_TRUE;
245674ae819SStefano Zampini       s++;
246674ae819SStefano Zampini     }
2472d456bbaSstefano_zampini     if (j < 1) skip_node = PETSC_TRUE;
2482d456bbaSstefano_zampini     if (mat_graph->special_dof[i] == PCBDDCGRAPH_DIRICHLET_MARK) skip_node = PETSC_TRUE;
249674ae819SStefano Zampini     if (!skip_node) {
250674ae819SStefano Zampini       if (fully_redundant) {
251674ae819SStefano Zampini         /* fully redundant set of lagrange multipliers */
252674ae819SStefano Zampini         n_lambda_for_dof = (j*(j+1))/2;
253674ae819SStefano Zampini       } else {
254674ae819SStefano Zampini         n_lambda_for_dof = j;
255674ae819SStefano Zampini       }
256674ae819SStefano Zampini       n_local_lambda += j;
257674ae819SStefano Zampini       /* needed to evaluate global number of lagrange multipliers */
258674ae819SStefano Zampini       array[i]=(1.0*n_lambda_for_dof)/(j+1.0); /* already scaled for the next global sum */
259674ae819SStefano Zampini       /* store some data needed */
260674ae819SStefano Zampini       dual_dofs_boundary_indices[partial_sum] = n_boundary_dofs-1;
261674ae819SStefano Zampini       aux_local_numbering_1[partial_sum] = i;
262674ae819SStefano Zampini       aux_local_numbering_2[partial_sum] = n_lambda_for_dof;
263674ae819SStefano Zampini       partial_sum++;
264674ae819SStefano Zampini     }
265674ae819SStefano Zampini   }
2669566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(pcis->vec1_N,&array));
2679566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isvert,&vertex_indices));
2689566063dSJacob Faibussowitsch   PetscCall(PCBDDCGraphRestoreCandidatesIS(pcbddc->mat_graph,NULL,NULL,NULL,NULL,&isvert));
2692d456bbaSstefano_zampini   dual_size = partial_sum;
270674ae819SStefano Zampini 
271674ae819SStefano Zampini   /* compute global ordering of lagrange multipliers and associate l2g map */
2729566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(comm,partial_sum,aux_local_numbering_1,PETSC_COPY_VALUES,&subset_n));
2739566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingApplyIS(pcis->mapping,subset_n,&subset));
2749566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&subset_n));
2759566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(comm,partial_sum,aux_local_numbering_2,PETSC_OWN_POINTER,&subset_mult));
2769566063dSJacob Faibussowitsch   PetscCall(ISRenumber(subset,subset_mult,&fetidpmat_ctx->n_lambda,&subset_n));
2779566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&subset));
2783d996552SStefano Zampini 
27976bd3646SJed Brown   if (PetscDefined(USE_DEBUG)) {
2809566063dSJacob Faibussowitsch     PetscCall(VecSet(pcis->vec1_global,0.0));
2819566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(matis->rctx,pcis->vec1_N,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE));
2829566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(matis->rctx,pcis->vec1_N,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE));
2839566063dSJacob Faibussowitsch     PetscCall(VecSum(pcis->vec1_global,&scalar_value));
2844fe826edSStefano Zampini     i = (PetscInt)PetscRealPart(scalar_value);
28563a3b9bcSJacob 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);
28676bd3646SJed Brown   }
287674ae819SStefano Zampini 
288674ae819SStefano Zampini   /* init data for scaling factors exchange */
28925d06dbeSstefano_zampini   if (!pcbddc->use_deluxe_scaling) {
29025d06dbeSstefano_zampini     PetscInt    *ptrs_buffer,neigh_position;
29125d06dbeSstefano_zampini     PetscScalar *send_buffer,*recv_buffer;
29225d06dbeSstefano_zampini     MPI_Request *send_reqs,*recv_reqs;
29325d06dbeSstefano_zampini 
294674ae819SStefano Zampini     partial_sum = 0;
2959566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(pcis->n_neigh,&ptrs_buffer));
2969566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(PetscMax(pcis->n_neigh-1,0),&send_reqs));
2979566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(PetscMax(pcis->n_neigh-1,0),&recv_reqs));
2989566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(pcis->n+1,&all_factors));
2994b2aedd3SStefano Zampini     if (pcis->n_neigh > 0) ptrs_buffer[0]=0;
300674ae819SStefano Zampini     for (i=1;i<pcis->n_neigh;i++) {
301674ae819SStefano Zampini       partial_sum += pcis->n_shared[i];
302674ae819SStefano Zampini       ptrs_buffer[i] = ptrs_buffer[i-1]+pcis->n_shared[i];
303674ae819SStefano Zampini     }
3049566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(partial_sum,&send_buffer));
3059566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(partial_sum,&recv_buffer));
3069566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(partial_sum,&all_factors[0]));
307674ae819SStefano Zampini     for (i=0;i<pcis->n-1;i++) {
308674ae819SStefano Zampini       j = mat_graph->count[i];
309674ae819SStefano Zampini       all_factors[i+1]=all_factors[i]+j;
310674ae819SStefano Zampini     }
31125d06dbeSstefano_zampini 
312674ae819SStefano Zampini     /* scatter B scaling to N vec */
3139566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(pcis->N_to_B,pcis->D,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE));
3149566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(pcis->N_to_B,pcis->D,pcis->vec1_N,INSERT_VALUES,SCATTER_REVERSE));
315674ae819SStefano Zampini     /* communications */
3169566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(pcis->vec1_N,(const PetscScalar**)&array));
317674ae819SStefano Zampini     for (i=1;i<pcis->n_neigh;i++) {
318674ae819SStefano Zampini       for (j=0;j<pcis->n_shared[i];j++) {
319674ae819SStefano Zampini         send_buffer[ptrs_buffer[i-1]+j]=array[pcis->shared[i][j]];
320674ae819SStefano Zampini       }
3219566063dSJacob Faibussowitsch       PetscCall(PetscMPIIntCast(ptrs_buffer[i]-ptrs_buffer[i-1],&buf_size));
3229566063dSJacob Faibussowitsch       PetscCall(PetscMPIIntCast(pcis->neigh[i],&neigh));
3239566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Isend(&send_buffer[ptrs_buffer[i-1]],buf_size,MPIU_SCALAR,neigh,0,comm,&send_reqs[i-1]));
3249566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Irecv(&recv_buffer[ptrs_buffer[i-1]],buf_size,MPIU_SCALAR,neigh,0,comm,&recv_reqs[i-1]));
325674ae819SStefano Zampini     }
3269566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(pcis->vec1_N,(const PetscScalar**)&array));
3274b2aedd3SStefano Zampini     if (pcis->n_neigh > 0) {
3289566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Waitall(pcis->n_neigh-1,recv_reqs,MPI_STATUSES_IGNORE));
3294b2aedd3SStefano Zampini     }
330674ae819SStefano Zampini     /* put values in correct places */
331674ae819SStefano Zampini     for (i=1;i<pcis->n_neigh;i++) {
332674ae819SStefano Zampini       for (j=0;j<pcis->n_shared[i];j++) {
333674ae819SStefano Zampini         k = pcis->shared[i][j];
334674ae819SStefano Zampini         neigh_position = 0;
335674ae819SStefano Zampini         while (mat_graph->neighbours_set[k][neigh_position] != pcis->neigh[i]) {neigh_position++;}
336674ae819SStefano Zampini         all_factors[k][neigh_position]=recv_buffer[ptrs_buffer[i-1]+j];
337674ae819SStefano Zampini       }
338674ae819SStefano Zampini     }
3394b2aedd3SStefano Zampini     if (pcis->n_neigh > 0) {
3409566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Waitall(pcis->n_neigh-1,send_reqs,MPI_STATUSES_IGNORE));
3414b2aedd3SStefano Zampini     }
3429566063dSJacob Faibussowitsch     PetscCall(PetscFree(send_reqs));
3439566063dSJacob Faibussowitsch     PetscCall(PetscFree(recv_reqs));
3449566063dSJacob Faibussowitsch     PetscCall(PetscFree(send_buffer));
3459566063dSJacob Faibussowitsch     PetscCall(PetscFree(recv_buffer));
3469566063dSJacob Faibussowitsch     PetscCall(PetscFree(ptrs_buffer));
34725d06dbeSstefano_zampini   }
348674ae819SStefano Zampini 
349674ae819SStefano Zampini   /* Compute B and B_delta (local actions) */
3509566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(pcis->n_neigh,&aux_sums));
3519566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n_local_lambda,&l2g_indices));
3529566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n_local_lambda,&vals_B_delta));
3539566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n_local_lambda,&cols_B_delta));
35425d06dbeSstefano_zampini   if (!pcbddc->use_deluxe_scaling) {
3559566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n_local_lambda,&scaling_factors));
35625d06dbeSstefano_zampini   } else {
35725d06dbeSstefano_zampini     scaling_factors = NULL;
35825d06dbeSstefano_zampini     all_factors     = NULL;
35925d06dbeSstefano_zampini   }
3609566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(subset_n,&aux_global_numbering));
361674ae819SStefano Zampini   partial_sum=0;
362dc456d91SStefano Zampini   cum = 0;
363674ae819SStefano Zampini   for (i=0;i<dual_size;i++) {
364dc456d91SStefano Zampini     n_global_lambda = aux_global_numbering[cum];
365674ae819SStefano Zampini     j = mat_graph->count[aux_local_numbering_1[i]];
366674ae819SStefano Zampini     aux_sums[0]=0;
367674ae819SStefano Zampini     for (s=1;s<j;s++) {
368674ae819SStefano Zampini       aux_sums[s]=aux_sums[s-1]+j-s+1;
369674ae819SStefano Zampini     }
37025d06dbeSstefano_zampini     if (all_factors) array = all_factors[aux_local_numbering_1[i]];
371674ae819SStefano Zampini     n_neg_values = 0;
3722a7da448SStefano Zampini     while (n_neg_values < j && mat_graph->neighbours_set[aux_local_numbering_1[i]][n_neg_values] < rank) {n_neg_values++;}
373674ae819SStefano Zampini     n_pos_values = j - n_neg_values;
374674ae819SStefano Zampini     if (fully_redundant) {
375674ae819SStefano Zampini       for (s=0;s<n_neg_values;s++) {
376674ae819SStefano Zampini         l2g_indices    [partial_sum+s]=aux_sums[s]+n_neg_values-s-1+n_global_lambda;
377674ae819SStefano Zampini         cols_B_delta   [partial_sum+s]=dual_dofs_boundary_indices[i];
378674ae819SStefano Zampini         vals_B_delta   [partial_sum+s]=-1.0;
37925d06dbeSstefano_zampini         if (!pcbddc->use_deluxe_scaling) scaling_factors[partial_sum+s]=array[s];
380674ae819SStefano Zampini       }
381674ae819SStefano Zampini       for (s=0;s<n_pos_values;s++) {
382674ae819SStefano Zampini         l2g_indices    [partial_sum+s+n_neg_values]=aux_sums[n_neg_values]+s+n_global_lambda;
383674ae819SStefano Zampini         cols_B_delta   [partial_sum+s+n_neg_values]=dual_dofs_boundary_indices[i];
384674ae819SStefano Zampini         vals_B_delta   [partial_sum+s+n_neg_values]=1.0;
38525d06dbeSstefano_zampini         if (!pcbddc->use_deluxe_scaling) scaling_factors[partial_sum+s+n_neg_values]=array[s+n_neg_values];
386674ae819SStefano Zampini       }
387674ae819SStefano Zampini       partial_sum += j;
388674ae819SStefano Zampini     } else {
389674ae819SStefano Zampini       /* l2g_indices and default cols and vals of B_delta */
390674ae819SStefano Zampini       for (s=0;s<j;s++) {
391674ae819SStefano Zampini         l2g_indices    [partial_sum+s]=n_global_lambda+s;
392674ae819SStefano Zampini         cols_B_delta   [partial_sum+s]=dual_dofs_boundary_indices[i];
393674ae819SStefano Zampini         vals_B_delta   [partial_sum+s]=0.0;
394674ae819SStefano Zampini       }
395674ae819SStefano Zampini       /* B_delta */
396674ae819SStefano Zampini       if (n_neg_values > 0) { /* there's a rank next to me to the left */
397674ae819SStefano Zampini         vals_B_delta   [partial_sum+n_neg_values-1]=-1.0;
398674ae819SStefano Zampini       }
399674ae819SStefano Zampini       if (n_neg_values < j) { /* there's a rank next to me to the right */
400674ae819SStefano Zampini         vals_B_delta   [partial_sum+n_neg_values]=1.0;
401674ae819SStefano Zampini       }
402674ae819SStefano Zampini       /* scaling as in Klawonn-Widlund 1999 */
40325d06dbeSstefano_zampini       if (!pcbddc->use_deluxe_scaling) {
404674ae819SStefano Zampini         for (s=0;s<n_neg_values;s++) {
405674ae819SStefano Zampini           scalar_value = 0.0;
40625d06dbeSstefano_zampini           for (k=0;k<s+1;k++) scalar_value += array[k];
407674ae819SStefano Zampini           scaling_factors[partial_sum+s] = -scalar_value;
408674ae819SStefano Zampini         }
409674ae819SStefano Zampini         for (s=0;s<n_pos_values;s++) {
410674ae819SStefano Zampini           scalar_value = 0.0;
41125d06dbeSstefano_zampini           for (k=s+n_neg_values;k<j;k++) scalar_value += array[k];
412674ae819SStefano Zampini           scaling_factors[partial_sum+s+n_neg_values] = scalar_value;
413674ae819SStefano Zampini         }
41425d06dbeSstefano_zampini       }
415674ae819SStefano Zampini       partial_sum += j;
416674ae819SStefano Zampini     }
417dc456d91SStefano Zampini     cum += aux_local_numbering_2[i];
418674ae819SStefano Zampini   }
4199566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(subset_n,&aux_global_numbering));
4209566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&subset_mult));
4219566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&subset_n));
4229566063dSJacob Faibussowitsch   PetscCall(PetscFree(aux_sums));
4239566063dSJacob Faibussowitsch   PetscCall(PetscFree(aux_local_numbering_1));
4249566063dSJacob Faibussowitsch   PetscCall(PetscFree(dual_dofs_boundary_indices));
42525d06dbeSstefano_zampini   if (all_factors) {
4269566063dSJacob Faibussowitsch     PetscCall(PetscFree(all_factors[0]));
4279566063dSJacob Faibussowitsch     PetscCall(PetscFree(all_factors));
42825d06dbeSstefano_zampini   }
429674ae819SStefano Zampini 
430674ae819SStefano Zampini   /* Create local part of B_delta */
4319566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_SELF,&fetidpmat_ctx->B_delta));
4329566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(fetidpmat_ctx->B_delta,n_local_lambda,pcis->n_B,n_local_lambda,pcis->n_B));
4339566063dSJacob Faibussowitsch   PetscCall(MatSetType(fetidpmat_ctx->B_delta,MATSEQAIJ));
4349566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(fetidpmat_ctx->B_delta,1,NULL));
4359566063dSJacob Faibussowitsch   PetscCall(MatSetOption(fetidpmat_ctx->B_delta,MAT_IGNORE_ZERO_ENTRIES,PETSC_TRUE));
436674ae819SStefano Zampini   for (i=0;i<n_local_lambda;i++) {
4379566063dSJacob Faibussowitsch     PetscCall(MatSetValue(fetidpmat_ctx->B_delta,i,cols_B_delta[i],vals_B_delta[i],INSERT_VALUES));
438674ae819SStefano Zampini   }
4399566063dSJacob Faibussowitsch   PetscCall(PetscFree(vals_B_delta));
4409566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(fetidpmat_ctx->B_delta,MAT_FINAL_ASSEMBLY));
4419566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(fetidpmat_ctx->B_delta,MAT_FINAL_ASSEMBLY));
442674ae819SStefano Zampini 
44325d06dbeSstefano_zampini   BD1 = NULL;
44425d06dbeSstefano_zampini   BD2 = NULL;
445674ae819SStefano Zampini   if (fully_redundant) {
44628b400f6SJacob Faibussowitsch     PetscCheck(!pcbddc->use_deluxe_scaling,comm,PETSC_ERR_SUP,"Deluxe FETIDP with fully-redundant multipliers to be implemented");
4479566063dSJacob Faibussowitsch     PetscCall(MatCreate(PETSC_COMM_SELF,&ScalingMat));
4489566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(ScalingMat,n_local_lambda,n_local_lambda,n_local_lambda,n_local_lambda));
4499566063dSJacob Faibussowitsch     PetscCall(MatSetType(ScalingMat,MATSEQAIJ));
4509566063dSJacob Faibussowitsch     PetscCall(MatSeqAIJSetPreallocation(ScalingMat,1,NULL));
451674ae819SStefano Zampini     for (i=0;i<n_local_lambda;i++) {
4529566063dSJacob Faibussowitsch       PetscCall(MatSetValue(ScalingMat,i,i,scaling_factors[i],INSERT_VALUES));
453674ae819SStefano Zampini     }
4549566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(ScalingMat,MAT_FINAL_ASSEMBLY));
4559566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(ScalingMat,MAT_FINAL_ASSEMBLY));
4569566063dSJacob Faibussowitsch     PetscCall(MatMatMult(ScalingMat,fetidpmat_ctx->B_delta,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&fetidpmat_ctx->B_Ddelta));
4579566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&ScalingMat));
458674ae819SStefano Zampini   } else {
4599566063dSJacob Faibussowitsch     PetscCall(MatCreate(PETSC_COMM_SELF,&fetidpmat_ctx->B_Ddelta));
4609566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(fetidpmat_ctx->B_Ddelta,n_local_lambda,pcis->n_B,n_local_lambda,pcis->n_B));
461524221abSStefano Zampini     if (!pcbddc->use_deluxe_scaling || !pcbddc->sub_schurs) {
4629566063dSJacob Faibussowitsch       PetscCall(MatSetType(fetidpmat_ctx->B_Ddelta,MATSEQAIJ));
4639566063dSJacob Faibussowitsch       PetscCall(MatSeqAIJSetPreallocation(fetidpmat_ctx->B_Ddelta,1,NULL));
464674ae819SStefano Zampini       for (i=0;i<n_local_lambda;i++) {
4659566063dSJacob Faibussowitsch         PetscCall(MatSetValue(fetidpmat_ctx->B_Ddelta,i,cols_B_delta[i],scaling_factors[i],INSERT_VALUES));
466674ae819SStefano Zampini       }
4679566063dSJacob Faibussowitsch       PetscCall(MatAssemblyBegin(fetidpmat_ctx->B_Ddelta,MAT_FINAL_ASSEMBLY));
4689566063dSJacob Faibussowitsch       PetscCall(MatAssemblyEnd(fetidpmat_ctx->B_Ddelta,MAT_FINAL_ASSEMBLY));
46925d06dbeSstefano_zampini     } else {
47025d06dbeSstefano_zampini       /* scaling as in Klawonn-Widlund 1999 */
47125d06dbeSstefano_zampini       PCBDDCDeluxeScaling deluxe_ctx = pcbddc->deluxe_ctx;
47225d06dbeSstefano_zampini       PCBDDCSubSchurs     sub_schurs = pcbddc->sub_schurs;
47325d06dbeSstefano_zampini       Mat                 T;
47425d06dbeSstefano_zampini       PetscScalar         *W,lwork,*Bwork;
475451a39c7SStefano Zampini       const PetscInt      *idxs = NULL;
47625d06dbeSstefano_zampini       PetscInt            cum,mss,*nnz;
47725d06dbeSstefano_zampini       PetscBLASInt        *pivots,B_lwork,B_N,B_ierr;
47825d06dbeSstefano_zampini 
47928b400f6SJacob Faibussowitsch       PetscCheck(pcbddc->deluxe_singlemat,comm,PETSC_ERR_USER,"Cannot compute B_Ddelta! rerun with -pc_bddc_deluxe_singlemat");
48025d06dbeSstefano_zampini       mss  = 0;
4819566063dSJacob Faibussowitsch       PetscCall(PetscCalloc1(pcis->n_B,&nnz));
48225d06dbeSstefano_zampini       if (sub_schurs->is_Ej_all) {
4839566063dSJacob Faibussowitsch         PetscCall(ISGetIndices(sub_schurs->is_Ej_all,&idxs));
48425d06dbeSstefano_zampini         for (i=0,cum=0;i<sub_schurs->n_subs;i++) {
48525d06dbeSstefano_zampini           PetscInt subset_size;
48625d06dbeSstefano_zampini 
4879566063dSJacob Faibussowitsch           PetscCall(ISGetLocalSize(sub_schurs->is_subs[i],&subset_size));
48825d06dbeSstefano_zampini           for (j=0;j<subset_size;j++) nnz[idxs[j+cum]] = subset_size;
48925d06dbeSstefano_zampini           mss  = PetscMax(mss,subset_size);
49025d06dbeSstefano_zampini           cum += subset_size;
49125d06dbeSstefano_zampini         }
49225d06dbeSstefano_zampini       }
4939566063dSJacob Faibussowitsch       PetscCall(MatCreate(PETSC_COMM_SELF,&T));
4949566063dSJacob Faibussowitsch       PetscCall(MatSetSizes(T,pcis->n_B,pcis->n_B,pcis->n_B,pcis->n_B));
4959566063dSJacob Faibussowitsch       PetscCall(MatSetType(T,MATSEQAIJ));
4969566063dSJacob Faibussowitsch       PetscCall(MatSeqAIJSetPreallocation(T,0,nnz));
4979566063dSJacob Faibussowitsch       PetscCall(PetscFree(nnz));
49825d06dbeSstefano_zampini 
49925d06dbeSstefano_zampini       /* workspace allocation */
5007f00194aSstefano_zampini       B_lwork = 0;
5017f00194aSstefano_zampini       if (mss) {
502be8bcbd0Sstefano_zampini         PetscScalar dummy = 1;
503be8bcbd0Sstefano_zampini 
50425d06dbeSstefano_zampini         B_lwork = -1;
5059566063dSJacob Faibussowitsch         PetscCall(PetscBLASIntCast(mss,&B_N));
5069566063dSJacob Faibussowitsch         PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
507*792fecdfSBarry Smith         PetscCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,&dummy,&B_N,&B_N,&lwork,&B_lwork,&B_ierr));
5089566063dSJacob Faibussowitsch         PetscCall(PetscFPTrapPop());
50928b400f6SJacob Faibussowitsch         PetscCheck(!B_ierr,PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to GETRI Lapack routine %d",(int)B_ierr);
5109566063dSJacob Faibussowitsch         PetscCall(PetscBLASIntCast((PetscInt)PetscRealPart(lwork),&B_lwork));
5117f00194aSstefano_zampini       }
5129566063dSJacob Faibussowitsch       PetscCall(PetscMalloc3(mss*mss,&W,mss,&pivots,B_lwork,&Bwork));
51325d06dbeSstefano_zampini 
51425d06dbeSstefano_zampini       for (i=0,cum=0;i<sub_schurs->n_subs;i++) {
5151683a169SBarry Smith         const PetscScalar *M;
51625d06dbeSstefano_zampini         PetscInt          subset_size;
51725d06dbeSstefano_zampini 
5189566063dSJacob Faibussowitsch         PetscCall(ISGetLocalSize(sub_schurs->is_subs[i],&subset_size));
5199566063dSJacob Faibussowitsch         PetscCall(PetscBLASIntCast(subset_size,&B_N));
5209566063dSJacob Faibussowitsch         PetscCall(MatDenseGetArrayRead(deluxe_ctx->seq_mat[i],&M));
5219566063dSJacob Faibussowitsch         PetscCall(PetscArraycpy(W,M,subset_size*subset_size));
5229566063dSJacob Faibussowitsch         PetscCall(MatDenseRestoreArrayRead(deluxe_ctx->seq_mat[i],&M));
5239566063dSJacob Faibussowitsch         PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
524*792fecdfSBarry Smith         PetscCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,W,&B_N,pivots,&B_ierr));
52528b400f6SJacob Faibussowitsch         PetscCheck(!B_ierr,PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
526*792fecdfSBarry Smith         PetscCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,W,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
52728b400f6SJacob Faibussowitsch         PetscCheck(!B_ierr,PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
5289566063dSJacob Faibussowitsch         PetscCall(PetscFPTrapPop());
529451a39c7SStefano Zampini         /* silent static analyzer */
53028b400f6SJacob Faibussowitsch         PetscCheck(idxs,PETSC_COMM_SELF,PETSC_ERR_PLIB,"IDXS not present");
5319566063dSJacob Faibussowitsch         PetscCall(MatSetValues(T,subset_size,idxs+cum,subset_size,idxs+cum,W,INSERT_VALUES));
53225d06dbeSstefano_zampini         cum += subset_size;
53325d06dbeSstefano_zampini       }
5349566063dSJacob Faibussowitsch       PetscCall(MatAssemblyBegin(T,MAT_FINAL_ASSEMBLY));
5359566063dSJacob Faibussowitsch       PetscCall(MatAssemblyEnd(T,MAT_FINAL_ASSEMBLY));
5369566063dSJacob Faibussowitsch       PetscCall(MatMatTransposeMult(T,fetidpmat_ctx->B_delta,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&BD1));
5379566063dSJacob Faibussowitsch       PetscCall(MatMatMult(fetidpmat_ctx->B_delta,BD1,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&BD2));
5389566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&T));
5399566063dSJacob Faibussowitsch       PetscCall(PetscFree3(W,pivots,Bwork));
54025d06dbeSstefano_zampini       if (sub_schurs->is_Ej_all) {
5419566063dSJacob Faibussowitsch         PetscCall(ISRestoreIndices(sub_schurs->is_Ej_all,&idxs));
54225d06dbeSstefano_zampini       }
54325d06dbeSstefano_zampini     }
544674ae819SStefano Zampini   }
5459566063dSJacob Faibussowitsch   PetscCall(PetscFree(scaling_factors));
5469566063dSJacob Faibussowitsch   PetscCall(PetscFree(cols_B_delta));
547674ae819SStefano Zampini 
548457d4a33Sstefano_zampini   /* Layout of multipliers */
5499566063dSJacob Faibussowitsch   PetscCall(PetscLayoutCreate(comm,&llay));
5509566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetBlockSize(llay,1));
5519566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetSize(llay,fetidpmat_ctx->n_lambda));
5529566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(llay));
5539566063dSJacob Faibussowitsch   PetscCall(PetscLayoutGetLocalSize(llay,&fetidpmat_ctx->n));
554457d4a33Sstefano_zampini 
555457d4a33Sstefano_zampini   /* Local work vector of multipliers */
5569566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_SELF,&fetidpmat_ctx->lambda_local));
5579566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(fetidpmat_ctx->lambda_local,n_local_lambda,n_local_lambda));
5589566063dSJacob Faibussowitsch   PetscCall(VecSetType(fetidpmat_ctx->lambda_local,VECSEQ));
559457d4a33Sstefano_zampini 
56025d06dbeSstefano_zampini   if (BD2) {
56125d06dbeSstefano_zampini     ISLocalToGlobalMapping l2g;
56225d06dbeSstefano_zampini     Mat                    T,TA,*pT;
56325d06dbeSstefano_zampini     IS                     is;
56425d06dbeSstefano_zampini     PetscInt               nl,N;
56525d06dbeSstefano_zampini     BDdelta_DN             ctx;
56625d06dbeSstefano_zampini 
5679566063dSJacob Faibussowitsch     PetscCall(PetscLayoutGetLocalSize(llay,&nl));
5689566063dSJacob Faibussowitsch     PetscCall(PetscLayoutGetSize(llay,&N));
5699566063dSJacob Faibussowitsch     PetscCall(MatCreate(comm,&T));
5709566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(T,nl,nl,N,N));
5719566063dSJacob Faibussowitsch     PetscCall(MatSetType(T,MATIS));
5729566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingCreate(comm,1,n_local_lambda,l2g_indices,PETSC_COPY_VALUES,&l2g));
5739566063dSJacob Faibussowitsch     PetscCall(MatSetLocalToGlobalMapping(T,l2g,l2g));
5749566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingDestroy(&l2g));
5759566063dSJacob Faibussowitsch     PetscCall(MatISSetLocalMat(T,BD2));
5769566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(T,MAT_FINAL_ASSEMBLY));
5779566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(T,MAT_FINAL_ASSEMBLY));
5789566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&BD2));
5799566063dSJacob Faibussowitsch     PetscCall(MatConvert(T,MATAIJ,MAT_INITIAL_MATRIX,&TA));
5809566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&T));
5819566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(comm,n_local_lambda,l2g_indices,PETSC_USE_POINTER,&is));
5829566063dSJacob Faibussowitsch     PetscCall(MatCreateSubMatrices(TA,1,&is,&is,MAT_INITIAL_MATRIX,&pT));
5839566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&TA));
5849566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&is));
58525d06dbeSstefano_zampini     BD2  = pT[0];
5869566063dSJacob Faibussowitsch     PetscCall(PetscFree(pT));
58725d06dbeSstefano_zampini 
58825d06dbeSstefano_zampini     /* B_Ddelta for non-redundant multipliers with deluxe scaling */
5899566063dSJacob Faibussowitsch     PetscCall(PetscNew(&ctx));
5909566063dSJacob Faibussowitsch     PetscCall(MatSetType(fetidpmat_ctx->B_Ddelta,MATSHELL));
5919566063dSJacob Faibussowitsch     PetscCall(MatShellSetContext(fetidpmat_ctx->B_Ddelta,ctx));
5929566063dSJacob Faibussowitsch     PetscCall(MatShellSetOperation(fetidpmat_ctx->B_Ddelta,MATOP_MULT,(void (*)(void))MatMult_BDdelta_deluxe_nonred));
5939566063dSJacob Faibussowitsch     PetscCall(MatShellSetOperation(fetidpmat_ctx->B_Ddelta,MATOP_MULT_TRANSPOSE,(void (*)(void))MatMultTranspose_BDdelta_deluxe_nonred));
5949566063dSJacob Faibussowitsch     PetscCall(MatShellSetOperation(fetidpmat_ctx->B_Ddelta,MATOP_DESTROY,(void (*)(void))MatDestroy_BDdelta_deluxe_nonred));
5959566063dSJacob Faibussowitsch     PetscCall(MatSetUp(fetidpmat_ctx->B_Ddelta));
59625d06dbeSstefano_zampini 
5979566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)BD1));
59825d06dbeSstefano_zampini     ctx->BD = BD1;
5999566063dSJacob Faibussowitsch     PetscCall(KSPCreate(PETSC_COMM_SELF,&ctx->kBD));
6009566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(ctx->kBD,BD2,BD2));
6019566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(fetidpmat_ctx->lambda_local,&ctx->work));
60225d06dbeSstefano_zampini     fetidpmat_ctx->deluxe_nonred = PETSC_TRUE;
60325d06dbeSstefano_zampini   }
6049566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&BD1));
6059566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&BD2));
60625d06dbeSstefano_zampini 
60725d06dbeSstefano_zampini   /* fetidpmat sizes */
60825d06dbeSstefano_zampini   fetidpmat_ctx->n += nPgl;
60925d06dbeSstefano_zampini   fetidpmat_ctx->N  = fetidpmat_ctx->n_lambda+nPg;
61025d06dbeSstefano_zampini 
611457d4a33Sstefano_zampini   /* Global vector for FETI-DP linear system */
6129566063dSJacob Faibussowitsch   PetscCall(VecCreate(comm,&fetidp_global));
6139566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(fetidp_global,fetidpmat_ctx->n,fetidpmat_ctx->N));
6149566063dSJacob Faibussowitsch   PetscCall(VecSetType(fetidp_global,VECMPI));
6159566063dSJacob Faibussowitsch   PetscCall(VecSetUp(fetidp_global));
616457d4a33Sstefano_zampini 
6179eec4de8Sstefano_zampini   /* Decide layout for fetidp dofs: if it is a saddle point problem
6189eec4de8Sstefano_zampini      pressure is ordered first in the local part of the global vector
619457d4a33Sstefano_zampini      of the FETI-DP linear system */
620457d4a33Sstefano_zampini   if (nPg) {
621e1214c54Sstefano_zampini     Vec            v;
622af140850Sstefano_zampini     IS             IS_l2g_p,ais;
623457d4a33Sstefano_zampini     PetscLayout    alay;
624457d4a33Sstefano_zampini     const PetscInt *idxs,*pranges,*aranges,*lranges;
625af140850Sstefano_zampini     PetscInt       *l2g_indices_p,rst;
626e1214c54Sstefano_zampini     PetscMPIInt    rank;
627457d4a33Sstefano_zampini 
6289566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nPl,&l2g_indices_p));
6299566063dSJacob Faibussowitsch     PetscCall(VecGetLayout(fetidp_global,&alay));
6309566063dSJacob Faibussowitsch     PetscCall(PetscLayoutGetRanges(alay,&aranges));
6319566063dSJacob Faibussowitsch     PetscCall(PetscLayoutGetRanges(play,&pranges));
6329566063dSJacob Faibussowitsch     PetscCall(PetscLayoutGetRanges(llay,&lranges));
633e1214c54Sstefano_zampini 
6349566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)fetidp_global),&rank));
6359566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PetscObjectComm((PetscObject)fetidp_global),pranges[rank+1]-pranges[rank],aranges[rank],1,&fetidpmat_ctx->pressure));
6369566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)fetidpmat_ctx->pressure,"F_P"));
6379566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PetscObjectComm((PetscObject)fetidp_global),lranges[rank+1]-lranges[rank],aranges[rank]+pranges[rank+1]-pranges[rank],1,&fetidpmat_ctx->lagrange));
6389566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)fetidpmat_ctx->lagrange,"F_L"));
6399566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingGetIndices(l2gmap_p,&idxs));
640af140850Sstefano_zampini     /* shift local to global indices for pressure */
641457d4a33Sstefano_zampini     for (i=0;i<nPl;i++) {
642131c27b5Sprj-       PetscMPIInt owner;
643457d4a33Sstefano_zampini 
6449566063dSJacob Faibussowitsch       PetscCall(PetscLayoutFindOwner(play,idxs[i],&owner));
645457d4a33Sstefano_zampini       l2g_indices_p[i] = idxs[i]-pranges[owner]+aranges[owner];
646457d4a33Sstefano_zampini     }
6479566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingRestoreIndices(l2gmap_p,&idxs));
6489566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(comm,nPl,l2g_indices_p,PETSC_OWN_POINTER,&IS_l2g_p));
649af140850Sstefano_zampini 
650e1214c54Sstefano_zampini     /* local to global scatter for pressure */
6519566063dSJacob Faibussowitsch     PetscCall(VecScatterCreate(fetidpmat_ctx->vP,NULL,fetidp_global,IS_l2g_p,&fetidpmat_ctx->l2g_p));
6529566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&IS_l2g_p));
653457d4a33Sstefano_zampini 
654e1214c54Sstefano_zampini     /* scatter for lagrange multipliers only */
6559566063dSJacob Faibussowitsch     PetscCall(VecCreate(comm,&v));
6569566063dSJacob Faibussowitsch     PetscCall(VecSetType(v,VECSTANDARD));
6579566063dSJacob Faibussowitsch     PetscCall(VecSetLayout(v,llay));
6589566063dSJacob Faibussowitsch     PetscCall(VecSetUp(v));
6599566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(comm,n_local_lambda,l2g_indices,PETSC_COPY_VALUES,&ais));
6609566063dSJacob Faibussowitsch     PetscCall(VecScatterCreate(fetidpmat_ctx->lambda_local,NULL,v,ais,&fetidpmat_ctx->l2g_lambda_only));
6619566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&ais));
6629566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&v));
663e1214c54Sstefano_zampini 
664af140850Sstefano_zampini     /* shift local to global indices for multipliers */
665457d4a33Sstefano_zampini     for (i=0;i<n_local_lambda;i++) {
666131c27b5Sprj-       PetscInt    ps;
667131c27b5Sprj-       PetscMPIInt owner;
668457d4a33Sstefano_zampini 
6699566063dSJacob Faibussowitsch       PetscCall(PetscLayoutFindOwner(llay,l2g_indices[i],&owner));
670457d4a33Sstefano_zampini       ps = pranges[owner+1]-pranges[owner];
671457d4a33Sstefano_zampini       l2g_indices[i] = l2g_indices[i]-lranges[owner]+aranges[owner]+ps;
672457d4a33Sstefano_zampini     }
673457d4a33Sstefano_zampini 
674e1214c54Sstefano_zampini     /* scatter from alldofs to pressures global fetidp vector */
6759566063dSJacob Faibussowitsch     PetscCall(PetscLayoutGetRange(alay,&rst,NULL));
6769566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(comm,nPgl,rst,1,&ais));
6779566063dSJacob Faibussowitsch     PetscCall(VecScatterCreate(pcis->vec1_global,pP,fetidp_global,ais,&fetidpmat_ctx->g2g_p));
6789566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&ais));
679457d4a33Sstefano_zampini   }
6809566063dSJacob Faibussowitsch   PetscCall(PetscLayoutDestroy(&llay));
6819566063dSJacob Faibussowitsch   PetscCall(PetscLayoutDestroy(&play));
6829566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(comm,n_local_lambda,l2g_indices,PETSC_OWN_POINTER,&IS_l2g_lambda));
683a1c0d0daSstefano_zampini 
6849cc7774eSstefano_zampini   /* scatter from local to global multipliers */
6859566063dSJacob Faibussowitsch   PetscCall(VecScatterCreate(fetidpmat_ctx->lambda_local,NULL,fetidp_global,IS_l2g_lambda,&fetidpmat_ctx->l2g_lambda));
6869566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&IS_l2g_lambda));
6879566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingDestroy(&l2gmap_p));
6889566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&fetidp_global));
689457d4a33Sstefano_zampini 
690a1c0d0daSstefano_zampini   /* Create some work vectors needed by fetidp */
6919566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(pcis->vec1_B,&fetidpmat_ctx->temp_solution_B));
6929566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(pcis->vec1_D,&fetidpmat_ctx->temp_solution_D));
693674ae819SStefano Zampini   PetscFunctionReturn(0);
694674ae819SStefano Zampini }
695674ae819SStefano Zampini 
696674ae819SStefano Zampini PetscErrorCode PCBDDCSetupFETIDPPCContext(Mat fetimat, FETIDPPC_ctx fetidppc_ctx)
697674ae819SStefano Zampini {
698674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
69925d06dbeSstefano_zampini   PC_BDDC        *pcbddc = (PC_BDDC*)fetidppc_ctx->pc->data;
70025d06dbeSstefano_zampini   PC_IS          *pcis = (PC_IS*)fetidppc_ctx->pc->data;
701f28b6018SStefano Zampini   PetscBool      lumped = PETSC_FALSE;
702674ae819SStefano Zampini 
703674ae819SStefano Zampini   PetscFunctionBegin;
7049566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(fetimat,&mat_ctx));
705674ae819SStefano Zampini   /* get references from objects created when setting up feti mat context */
7069566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)mat_ctx->lambda_local));
707674ae819SStefano Zampini   fetidppc_ctx->lambda_local = mat_ctx->lambda_local;
7089566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)mat_ctx->B_Ddelta));
709674ae819SStefano Zampini   fetidppc_ctx->B_Ddelta = mat_ctx->B_Ddelta;
71025d06dbeSstefano_zampini   if (mat_ctx->deluxe_nonred) {
71125d06dbeSstefano_zampini     PC               pc,mpc;
71225d06dbeSstefano_zampini     BDdelta_DN       ctx;
7133ca39a21SBarry Smith     MatSolverType    solver;
71425d06dbeSstefano_zampini     const char       *prefix;
71525d06dbeSstefano_zampini 
7169566063dSJacob Faibussowitsch     PetscCall(MatShellGetContext(mat_ctx->B_Ddelta,&ctx));
7179566063dSJacob Faibussowitsch     PetscCall(KSPSetType(ctx->kBD,KSPPREONLY));
7189566063dSJacob Faibussowitsch     PetscCall(KSPGetPC(ctx->kBD,&mpc));
7199566063dSJacob Faibussowitsch     PetscCall(KSPGetPC(pcbddc->ksp_D,&pc));
7209566063dSJacob Faibussowitsch     PetscCall(PCSetType(mpc,PCLU));
7219566063dSJacob Faibussowitsch     PetscCall(PCFactorGetMatSolverType(pc,(MatSolverType*)&solver));
7221baa6e33SBarry Smith     if (solver) PetscCall(PCFactorSetMatSolverType(mpc,solver));
7239566063dSJacob Faibussowitsch     PetscCall(MatGetOptionsPrefix(fetimat,&prefix));
7249566063dSJacob Faibussowitsch     PetscCall(KSPSetOptionsPrefix(ctx->kBD,prefix));
7259566063dSJacob Faibussowitsch     PetscCall(KSPAppendOptionsPrefix(ctx->kBD,"bddelta_"));
7269566063dSJacob Faibussowitsch     PetscCall(KSPSetFromOptions(ctx->kBD));
72725d06dbeSstefano_zampini   }
72825d06dbeSstefano_zampini 
729e1214c54Sstefano_zampini   if (mat_ctx->l2g_lambda_only) {
7309566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)mat_ctx->l2g_lambda_only));
731e1214c54Sstefano_zampini     fetidppc_ctx->l2g_lambda = mat_ctx->l2g_lambda_only;
732e1214c54Sstefano_zampini   } else {
7339566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)mat_ctx->l2g_lambda));
734674ae819SStefano Zampini     fetidppc_ctx->l2g_lambda = mat_ctx->l2g_lambda;
735e1214c54Sstefano_zampini   }
736f28b6018SStefano Zampini   /* Dirichlet preconditioner */
7379566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL,((PetscObject)fetimat)->prefix,"-pc_lumped",&lumped,NULL));
738f28b6018SStefano Zampini   if (!lumped) {
7399feb908dSstefano_zampini     IS        iV;
7409c2d02cdSstefano_zampini     PetscBool discrete_harmonic = PETSC_FALSE;
7419c2d02cdSstefano_zampini 
7429566063dSJacob Faibussowitsch     PetscCall(PetscObjectQuery((PetscObject)fetidppc_ctx->pc,"__KSPFETIDP_iV",(PetscObject*)&iV));
7439feb908dSstefano_zampini     if (iV) {
7449566063dSJacob Faibussowitsch       PetscCall(PetscOptionsGetBool(NULL,((PetscObject)fetimat)->prefix,"-pc_discrete_harmonic",&discrete_harmonic,NULL));
7459c2d02cdSstefano_zampini     }
7469c2d02cdSstefano_zampini     if (discrete_harmonic) {
7479c2d02cdSstefano_zampini       KSP             sksp;
7489c2d02cdSstefano_zampini       PC              pc;
749111315fdSstefano_zampini       PCBDDCSubSchurs sub_schurs = pcbddc->sub_schurs;
7509c2d02cdSstefano_zampini       Mat             A_II,A_IB,A_BI;
751111315fdSstefano_zampini       IS              iP = NULL;
752111315fdSstefano_zampini       PetscBool       isshell,reuse = PETSC_FALSE;
7539c2d02cdSstefano_zampini       KSPType         ksptype;
7549c2d02cdSstefano_zampini       const char      *prefix;
7559c2d02cdSstefano_zampini 
7569c2d02cdSstefano_zampini       /*
7579c2d02cdSstefano_zampini         We constructs a Schur complement for
7589c2d02cdSstefano_zampini 
7599c2d02cdSstefano_zampini         | A_II A_ID |
7609c2d02cdSstefano_zampini         | A_DI A_DD |
7619c2d02cdSstefano_zampini 
7629c2d02cdSstefano_zampini         instead of
7639c2d02cdSstefano_zampini 
7649c2d02cdSstefano_zampini         | A_II  B^t_II A_ID |
7659c2d02cdSstefano_zampini         | B_II -C_II   B_ID |
7669c2d02cdSstefano_zampini         | A_DI  B^t_ID A_DD |
7679c2d02cdSstefano_zampini 
7689c2d02cdSstefano_zampini       */
769111315fdSstefano_zampini       if (sub_schurs && sub_schurs->reuse_solver) {
7709566063dSJacob Faibussowitsch         PetscCall(PetscObjectQuery((PetscObject)sub_schurs->A,"__KSPFETIDP_iP",(PetscObject*)&iP));
771111315fdSstefano_zampini         if (iP) reuse = PETSC_TRUE;
772111315fdSstefano_zampini       }
773111315fdSstefano_zampini       if (!reuse) {
774111315fdSstefano_zampini         IS       aB;
775111315fdSstefano_zampini         PetscInt nb;
7769566063dSJacob Faibussowitsch         PetscCall(ISGetLocalSize(pcis->is_B_local,&nb));
7779566063dSJacob Faibussowitsch         PetscCall(ISCreateStride(PetscObjectComm((PetscObject)pcis->A_II),nb,0,1,&aB));
7789566063dSJacob Faibussowitsch         PetscCall(MatCreateSubMatrix(pcis->A_II,iV,iV,MAT_INITIAL_MATRIX,&A_II));
7799566063dSJacob Faibussowitsch         PetscCall(MatCreateSubMatrix(pcis->A_IB,iV,aB,MAT_INITIAL_MATRIX,&A_IB));
7809566063dSJacob Faibussowitsch         PetscCall(MatCreateSubMatrix(pcis->A_BI,aB,iV,MAT_INITIAL_MATRIX,&A_BI));
7819566063dSJacob Faibussowitsch         PetscCall(ISDestroy(&aB));
782111315fdSstefano_zampini       } else {
7839566063dSJacob Faibussowitsch         PetscCall(MatCreateSubMatrix(sub_schurs->A,pcis->is_I_local,pcis->is_B_local,MAT_INITIAL_MATRIX,&A_IB));
7849566063dSJacob Faibussowitsch         PetscCall(MatCreateSubMatrix(sub_schurs->A,pcis->is_B_local,pcis->is_I_local,MAT_INITIAL_MATRIX,&A_BI));
7859566063dSJacob Faibussowitsch         PetscCall(PetscObjectReference((PetscObject)pcis->A_II));
786111315fdSstefano_zampini         A_II = pcis->A_II;
787111315fdSstefano_zampini       }
7889566063dSJacob Faibussowitsch       PetscCall(MatCreateSchurComplement(A_II,A_II,A_IB,A_BI,pcis->A_BB,&fetidppc_ctx->S_j));
7899c2d02cdSstefano_zampini 
7909c2d02cdSstefano_zampini       /* propagate settings of solver */
7919566063dSJacob Faibussowitsch       PetscCall(MatSchurComplementGetKSP(fetidppc_ctx->S_j,&sksp));
7929566063dSJacob Faibussowitsch       PetscCall(KSPGetType(pcis->ksp_D,&ksptype));
7939566063dSJacob Faibussowitsch       PetscCall(KSPSetType(sksp,ksptype));
7949566063dSJacob Faibussowitsch       PetscCall(KSPGetPC(pcis->ksp_D,&pc));
7959566063dSJacob Faibussowitsch       PetscCall(PetscObjectTypeCompare((PetscObject)pc,PCSHELL,&isshell));
7965334bea6Sstefano_zampini       if (!isshell) {
7973ca39a21SBarry Smith         MatSolverType    solver;
7985334bea6Sstefano_zampini         PCType           pctype;
7995334bea6Sstefano_zampini 
8009566063dSJacob Faibussowitsch         PetscCall(PCGetType(pc,&pctype));
8019566063dSJacob Faibussowitsch         PetscCall(PCFactorGetMatSolverType(pc,(MatSolverType*)&solver));
8029566063dSJacob Faibussowitsch         PetscCall(KSPGetPC(sksp,&pc));
8039566063dSJacob Faibussowitsch         PetscCall(PCSetType(pc,pctype));
8041baa6e33SBarry Smith         if (solver) PetscCall(PCFactorSetMatSolverType(pc,solver));
8055334bea6Sstefano_zampini       } else {
8069566063dSJacob Faibussowitsch         PetscCall(KSPGetPC(sksp,&pc));
8079566063dSJacob Faibussowitsch         PetscCall(PCSetType(pc,PCLU));
8085334bea6Sstefano_zampini       }
8099566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&A_II));
8109566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&A_IB));
8119566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&A_BI));
8129566063dSJacob Faibussowitsch       PetscCall(MatGetOptionsPrefix(fetimat,&prefix));
8139566063dSJacob Faibussowitsch       PetscCall(KSPSetOptionsPrefix(sksp,prefix));
8149566063dSJacob Faibussowitsch       PetscCall(KSPAppendOptionsPrefix(sksp,"harmonic_"));
8159566063dSJacob Faibussowitsch       PetscCall(KSPSetFromOptions(sksp));
816111315fdSstefano_zampini       if (reuse) {
8179566063dSJacob Faibussowitsch         PetscCall(KSPSetPC(sksp,sub_schurs->reuse_solver->interior_solver));
8189566063dSJacob Faibussowitsch         PetscCall(PetscObjectIncrementTabLevel((PetscObject)sub_schurs->reuse_solver->interior_solver,(PetscObject)sksp,0));
819111315fdSstefano_zampini       }
8209c2d02cdSstefano_zampini     } else { /* default Dirichlet preconditioner is pde-harmonic */
8219566063dSJacob Faibussowitsch       PetscCall(MatCreateSchurComplement(pcis->A_II,pcis->A_II,pcis->A_IB,pcis->A_BI,pcis->A_BB,&fetidppc_ctx->S_j));
8229566063dSJacob Faibussowitsch       PetscCall(MatSchurComplementSetKSP(fetidppc_ctx->S_j,pcis->ksp_D));
8239c2d02cdSstefano_zampini     }
824f28b6018SStefano Zampini   } else {
8259566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)pcis->A_BB));
826f28b6018SStefano Zampini     fetidppc_ctx->S_j = pcis->A_BB;
827f28b6018SStefano Zampini   }
828af140850Sstefano_zampini   /* saddle-point */
829af140850Sstefano_zampini   if (mat_ctx->xPg) {
8309566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)mat_ctx->xPg));
831af140850Sstefano_zampini     fetidppc_ctx->xPg = mat_ctx->xPg;
8329566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)mat_ctx->yPg));
833af140850Sstefano_zampini     fetidppc_ctx->yPg = mat_ctx->yPg;
834af140850Sstefano_zampini   }
835674ae819SStefano Zampini   PetscFunctionReturn(0);
836674ae819SStefano Zampini }
837674ae819SStefano Zampini 
838ef425aeaSstefano_zampini PetscErrorCode FETIDPMatMult_Kernel(Mat fetimat, Vec x, Vec y, PetscBool trans)
839674ae819SStefano Zampini {
840674ae819SStefano Zampini   FETIDPMat_ctx  mat_ctx;
841617d11aeSStefano Zampini   PC_BDDC        *pcbddc;
842674ae819SStefano Zampini   PC_IS          *pcis;
843674ae819SStefano Zampini 
844674ae819SStefano Zampini   PetscFunctionBegin;
8459566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(fetimat,&mat_ctx));
846674ae819SStefano Zampini   pcis = (PC_IS*)mat_ctx->pc->data;
847617d11aeSStefano Zampini   pcbddc = (PC_BDDC*)mat_ctx->pc->data;
848674ae819SStefano Zampini   /* Application of B_delta^T */
8499566063dSJacob Faibussowitsch   PetscCall(VecSet(pcis->vec1_B,0.));
8509566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(mat_ctx->l2g_lambda,x,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE));
8519566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(mat_ctx->l2g_lambda,x,mat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE));
8529566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(mat_ctx->B_delta,mat_ctx->lambda_local,pcis->vec1_B));
853af140850Sstefano_zampini 
854af140850Sstefano_zampini   /* Add contribution from saddle point */
855af140850Sstefano_zampini   if (mat_ctx->l2g_p) {
8569566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(mat_ctx->l2g_p,x,mat_ctx->vP,INSERT_VALUES,SCATTER_REVERSE));
8579566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(mat_ctx->l2g_p,x,mat_ctx->vP,INSERT_VALUES,SCATTER_REVERSE));
858af140850Sstefano_zampini     if (pcbddc->switch_static) {
859ef425aeaSstefano_zampini       if (trans) {
8609566063dSJacob Faibussowitsch         PetscCall(MatMultTranspose(mat_ctx->B_BI,mat_ctx->vP,pcis->vec1_D));
861ef425aeaSstefano_zampini       } else {
8629566063dSJacob Faibussowitsch         PetscCall(MatMult(mat_ctx->Bt_BI,mat_ctx->vP,pcis->vec1_D));
863af140850Sstefano_zampini       }
864ef425aeaSstefano_zampini     }
865ef425aeaSstefano_zampini     if (trans) {
8669566063dSJacob Faibussowitsch       PetscCall(MatMultTransposeAdd(mat_ctx->B_BB,mat_ctx->vP,pcis->vec1_B,pcis->vec1_B));
867ef425aeaSstefano_zampini     } else {
8689566063dSJacob Faibussowitsch       PetscCall(MatMultAdd(mat_ctx->Bt_BB,mat_ctx->vP,pcis->vec1_B,pcis->vec1_B));
869ef425aeaSstefano_zampini     }
870af140850Sstefano_zampini   } else {
8711baa6e33SBarry Smith     if (pcbddc->switch_static) PetscCall(VecSet(pcis->vec1_D,0.0));
872af140850Sstefano_zampini   }
873af140850Sstefano_zampini   /* Application of \widetilde{S}^-1 */
8749566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(pcbddc->benign_p0,pcbddc->benign_n));
8759566063dSJacob Faibussowitsch   PetscCall(PCBDDCApplyInterfacePreconditioner(mat_ctx->pc,trans));
8769566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(pcbddc->benign_p0,pcbddc->benign_n));
8779566063dSJacob Faibussowitsch   PetscCall(VecSet(y,0.0));
878674ae819SStefano Zampini   /* Application of B_delta */
8799566063dSJacob Faibussowitsch   PetscCall(MatMult(mat_ctx->B_delta,pcis->vec1_B,mat_ctx->lambda_local));
880af140850Sstefano_zampini   /* Contribution from boundary pressures */
881af140850Sstefano_zampini   if (mat_ctx->C) {
882af140850Sstefano_zampini     const PetscScalar *lx;
883af140850Sstefano_zampini     PetscScalar       *ly;
884af140850Sstefano_zampini 
88515579a77SStefano Zampini     /* pressure ordered first in the local part of x and y */
8869566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(x,&lx));
8879566063dSJacob Faibussowitsch     PetscCall(VecGetArray(y,&ly));
8889566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(mat_ctx->xPg,lx));
8899566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(mat_ctx->yPg,ly));
890ef425aeaSstefano_zampini     if (trans) {
8919566063dSJacob Faibussowitsch       PetscCall(MatMultTranspose(mat_ctx->C,mat_ctx->xPg,mat_ctx->yPg));
892ef425aeaSstefano_zampini     } else {
8939566063dSJacob Faibussowitsch       PetscCall(MatMult(mat_ctx->C,mat_ctx->xPg,mat_ctx->yPg));
894ef425aeaSstefano_zampini     }
8959566063dSJacob Faibussowitsch     PetscCall(VecResetArray(mat_ctx->xPg));
8969566063dSJacob Faibussowitsch     PetscCall(VecResetArray(mat_ctx->yPg));
8979566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(x,&lx));
8989566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(y,&ly));
899af140850Sstefano_zampini   }
900af140850Sstefano_zampini   /* Add contribution from saddle point */
901af140850Sstefano_zampini   if (mat_ctx->l2g_p) {
902ef425aeaSstefano_zampini     if (trans) {
9039566063dSJacob Faibussowitsch       PetscCall(MatMultTranspose(mat_ctx->Bt_BB,pcis->vec1_B,mat_ctx->vP));
904ef425aeaSstefano_zampini     } else {
9059566063dSJacob Faibussowitsch       PetscCall(MatMult(mat_ctx->B_BB,pcis->vec1_B,mat_ctx->vP));
906ef425aeaSstefano_zampini     }
907af140850Sstefano_zampini     if (pcbddc->switch_static) {
908ef425aeaSstefano_zampini       if (trans) {
9099566063dSJacob Faibussowitsch         PetscCall(MatMultTransposeAdd(mat_ctx->Bt_BI,pcis->vec1_D,mat_ctx->vP,mat_ctx->vP));
910ef425aeaSstefano_zampini       } else {
9119566063dSJacob Faibussowitsch         PetscCall(MatMultAdd(mat_ctx->B_BI,pcis->vec1_D,mat_ctx->vP,mat_ctx->vP));
912af140850Sstefano_zampini       }
913ef425aeaSstefano_zampini     }
9149566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(mat_ctx->l2g_p,mat_ctx->vP,y,ADD_VALUES,SCATTER_FORWARD));
9159566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(mat_ctx->l2g_p,mat_ctx->vP,y,ADD_VALUES,SCATTER_FORWARD));
916af140850Sstefano_zampini   }
9179566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(mat_ctx->l2g_lambda,mat_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD));
9189566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(mat_ctx->l2g_lambda,mat_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD));
919674ae819SStefano Zampini   PetscFunctionReturn(0);
920674ae819SStefano Zampini }
921674ae819SStefano Zampini 
922ef425aeaSstefano_zampini PetscErrorCode FETIDPMatMult(Mat fetimat, Vec x, Vec y)
923edf7251bSStefano Zampini {
924edf7251bSStefano Zampini   PetscFunctionBegin;
9259566063dSJacob Faibussowitsch   PetscCall(FETIDPMatMult_Kernel(fetimat,x,y,PETSC_FALSE));
926ef425aeaSstefano_zampini   PetscFunctionReturn(0);
927ef425aeaSstefano_zampini }
928ef425aeaSstefano_zampini 
929ef425aeaSstefano_zampini PetscErrorCode FETIDPMatMultTranspose(Mat fetimat, Vec x, Vec y)
930ef425aeaSstefano_zampini {
931ef425aeaSstefano_zampini   PetscFunctionBegin;
9329566063dSJacob Faibussowitsch   PetscCall(FETIDPMatMult_Kernel(fetimat,x,y,PETSC_TRUE));
933edf7251bSStefano Zampini   PetscFunctionReturn(0);
934edf7251bSStefano Zampini }
935edf7251bSStefano Zampini 
93616597391Sstefano_zampini PetscErrorCode FETIDPPCApply_Kernel(PC fetipc, Vec x, Vec y, PetscBool trans)
937674ae819SStefano Zampini {
938674ae819SStefano Zampini   FETIDPPC_ctx   pc_ctx;
939674ae819SStefano Zampini   PC_IS          *pcis;
940674ae819SStefano Zampini 
941674ae819SStefano Zampini   PetscFunctionBegin;
9429566063dSJacob Faibussowitsch   PetscCall(PCShellGetContext(fetipc,&pc_ctx));
943674ae819SStefano Zampini   pcis = (PC_IS*)pc_ctx->pc->data;
944674ae819SStefano Zampini   /* Application of B_Ddelta^T */
9459566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(pc_ctx->l2g_lambda,x,pc_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE));
9469566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(pc_ctx->l2g_lambda,x,pc_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE));
9479566063dSJacob Faibussowitsch   PetscCall(VecSet(pcis->vec2_B,0.0));
9489566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(pc_ctx->B_Ddelta,pc_ctx->lambda_local,pcis->vec2_B));
949ed6c3d69SStefano Zampini   /* Application of local Schur complement */
95016597391Sstefano_zampini   if (trans) {
9519566063dSJacob Faibussowitsch     PetscCall(MatMultTranspose(pc_ctx->S_j,pcis->vec2_B,pcis->vec1_B));
95216597391Sstefano_zampini   } else {
9539566063dSJacob Faibussowitsch     PetscCall(MatMult(pc_ctx->S_j,pcis->vec2_B,pcis->vec1_B));
95416597391Sstefano_zampini   }
955edf7251bSStefano Zampini   /* Application of B_Ddelta */
9569566063dSJacob Faibussowitsch   PetscCall(MatMult(pc_ctx->B_Ddelta,pcis->vec1_B,pc_ctx->lambda_local));
9579566063dSJacob Faibussowitsch   PetscCall(VecSet(y,0.0));
9589566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(pc_ctx->l2g_lambda,pc_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD));
9599566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(pc_ctx->l2g_lambda,pc_ctx->lambda_local,y,ADD_VALUES,SCATTER_FORWARD));
960edf7251bSStefano Zampini   PetscFunctionReturn(0);
961edf7251bSStefano Zampini }
962edf7251bSStefano Zampini 
96316597391Sstefano_zampini PetscErrorCode FETIDPPCApply(PC pc, Vec x, Vec y)
964edf7251bSStefano Zampini {
965edf7251bSStefano Zampini   PetscFunctionBegin;
9669566063dSJacob Faibussowitsch   PetscCall(FETIDPPCApply_Kernel(pc,x,y,PETSC_FALSE));
96716597391Sstefano_zampini   PetscFunctionReturn(0);
96816597391Sstefano_zampini }
96916597391Sstefano_zampini 
97016597391Sstefano_zampini PetscErrorCode FETIDPPCApplyTranspose(PC pc, Vec x, Vec y)
97116597391Sstefano_zampini {
97216597391Sstefano_zampini   PetscFunctionBegin;
9739566063dSJacob Faibussowitsch   PetscCall(FETIDPPCApply_Kernel(pc,x,y,PETSC_TRUE));
974674ae819SStefano Zampini   PetscFunctionReturn(0);
975674ae819SStefano Zampini }
976c45b8d2dSstefano_zampini 
977c45b8d2dSstefano_zampini PetscErrorCode FETIDPPCView(PC pc, PetscViewer viewer)
978c45b8d2dSstefano_zampini {
979c45b8d2dSstefano_zampini   FETIDPPC_ctx      pc_ctx;
980c45b8d2dSstefano_zampini   PetscBool         iascii;
981c45b8d2dSstefano_zampini   PetscViewer       sviewer;
982c45b8d2dSstefano_zampini 
983c45b8d2dSstefano_zampini   PetscFunctionBegin;
9849566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii));
985c45b8d2dSstefano_zampini   if (iascii) {
986c45b8d2dSstefano_zampini     PetscMPIInt rank;
98725d06dbeSstefano_zampini     PetscBool   isschur,isshell;
988c45b8d2dSstefano_zampini 
9899566063dSJacob Faibussowitsch     PetscCall(PCShellGetContext(pc,&pc_ctx));
9909566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank));
9919566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)pc_ctx->S_j,MATSCHURCOMPLEMENT,&isschur));
992c45b8d2dSstefano_zampini     if (isschur) {
9939566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"  Dirichlet preconditioner (just from rank 0)\n"));
994c45b8d2dSstefano_zampini     } else {
9959566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"  Lumped preconditioner (just from rank 0)\n"));
996c45b8d2dSstefano_zampini     }
9979566063dSJacob Faibussowitsch     PetscCall(PetscViewerGetSubViewer(viewer,PetscObjectComm((PetscObject)pc_ctx->S_j),&sviewer));
998dd400576SPatrick Sanan     if (rank == 0) {
9999566063dSJacob Faibussowitsch       PetscCall(PetscViewerPushFormat(sviewer,PETSC_VIEWER_ASCII_INFO));
10009566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPushTab(sviewer));
10019566063dSJacob Faibussowitsch       PetscCall(MatView(pc_ctx->S_j,sviewer));
10029566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPopTab(sviewer));
10039566063dSJacob Faibussowitsch       PetscCall(PetscViewerPopFormat(sviewer));
1004c45b8d2dSstefano_zampini     }
10059566063dSJacob Faibussowitsch     PetscCall(PetscViewerRestoreSubViewer(viewer,PetscObjectComm((PetscObject)pc_ctx->S_j),&sviewer));
10069566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)pc_ctx->B_Ddelta,MATSHELL,&isshell));
100725d06dbeSstefano_zampini     if (isshell) {
100825d06dbeSstefano_zampini       BDdelta_DN ctx;
10099566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"  FETI-DP BDdelta: DB^t * (B D^-1 B^t)^-1 for deluxe scaling (just from rank 0)\n"));
10109566063dSJacob Faibussowitsch       PetscCall(MatShellGetContext(pc_ctx->B_Ddelta,&ctx));
10119566063dSJacob Faibussowitsch       PetscCall(PetscViewerGetSubViewer(viewer,PetscObjectComm((PetscObject)pc_ctx->S_j),&sviewer));
1012dd400576SPatrick Sanan       if (rank == 0) {
101315579a77SStefano Zampini         PetscInt tl;
101415579a77SStefano Zampini 
10159566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIGetTab(sviewer,&tl));
10169566063dSJacob Faibussowitsch         PetscCall(PetscObjectSetTabLevel((PetscObject)ctx->kBD,tl));
10179566063dSJacob Faibussowitsch         PetscCall(KSPView(ctx->kBD,sviewer));
10189566063dSJacob Faibussowitsch         PetscCall(PetscViewerPushFormat(sviewer,PETSC_VIEWER_ASCII_INFO));
10199566063dSJacob Faibussowitsch         PetscCall(MatView(ctx->BD,sviewer));
10209566063dSJacob Faibussowitsch         PetscCall(PetscViewerPopFormat(sviewer));
102125d06dbeSstefano_zampini       }
10229566063dSJacob Faibussowitsch       PetscCall(PetscViewerRestoreSubViewer(viewer,PetscObjectComm((PetscObject)pc_ctx->S_j),&sviewer));
1023e5afcf28SBarry Smith     }
10249566063dSJacob Faibussowitsch     PetscCall(PetscViewerFlush(viewer));
1025c45b8d2dSstefano_zampini   }
1026c45b8d2dSstefano_zampini   PetscFunctionReturn(0);
1027c45b8d2dSstefano_zampini }
1028