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