xref: /petsc/src/ksp/pc/impls/bddc/bddcschurs.c (revision 4a6c6b0d6791c089c3142e43ded9f3ea81bca07c)
1 #include <../src/ksp/pc/impls/bddc/bddc.h>
2 #include <../src/ksp/pc/impls/bddc/bddcprivate.h>
3 #include <petscblaslapack.h>
4 
5 PETSC_STATIC_INLINE PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt*,PetscInt,PetscBT,PetscInt*,PetscInt*,PetscInt*);
6 static PetscErrorCode PCBDDCComputeExplicitSchur(Mat,PetscBool,MatReuse,Mat*);
7 static PetscErrorCode PCBDDCMumpsInteriorSolve(PC,Vec,Vec);
8 static PetscErrorCode PCBDDCMumpsCorrectionSolve(PC,Vec,Vec);
9 
10 #undef __FUNCT__
11 #define __FUNCT__ "PCBDDCMumpsCorrectionSolve_Private"
12 static PetscErrorCode PCBDDCMumpsCorrectionSolve_Private(PC pc, Vec rhs, Vec sol, PetscBool transpose)
13 {
14   PCBDDCReuseMumps ctx;
15 #if defined(PETSC_HAVE_MUMPS)
16   PetscInt         ival;
17 #endif
18   PetscErrorCode   ierr;
19 
20   PetscFunctionBegin;
21   ierr = PCShellGetContext(pc,(void **)&ctx);CHKERRQ(ierr);
22 #if defined(PETSC_HAVE_MUMPS)
23   ierr = MatMumpsGetIcntl(ctx->F,26,&ival);CHKERRQ(ierr);
24   ierr = MatMumpsSetIcntl(ctx->F,26,-1);CHKERRQ(ierr);
25 #endif
26   if (transpose) {
27     ierr = MatSolveTranspose(ctx->F,rhs,sol);CHKERRQ(ierr);
28   } else {
29     ierr = MatSolve(ctx->F,rhs,sol);CHKERRQ(ierr);
30   }
31 #if defined(PETSC_HAVE_MUMPS)
32   ierr = MatMumpsSetIcntl(ctx->F,26,ival);CHKERRQ(ierr);
33 #endif
34   PetscFunctionReturn(0);
35 }
36 
37 #undef __FUNCT__
38 #define __FUNCT__ "PCBDDCMumpsCorrectionSolve"
39 static PetscErrorCode PCBDDCMumpsCorrectionSolve(PC pc, Vec rhs, Vec sol)
40 {
41   PetscErrorCode   ierr;
42 
43   PetscFunctionBegin;
44   ierr = PCBDDCMumpsCorrectionSolve_Private(pc,rhs,sol,PETSC_FALSE);CHKERRQ(ierr);
45   PetscFunctionReturn(0);
46 }
47 
48 #undef __FUNCT__
49 #define __FUNCT__ "PCBDDCMumpsCorrectionSolveTranspose"
50 static PetscErrorCode PCBDDCMumpsCorrectionSolveTranspose(PC pc, Vec rhs, Vec sol)
51 {
52   PetscErrorCode   ierr;
53 
54   PetscFunctionBegin;
55   ierr = PCBDDCMumpsCorrectionSolve_Private(pc,rhs,sol,PETSC_TRUE);CHKERRQ(ierr);
56   PetscFunctionReturn(0);
57 }
58 
59 #undef __FUNCT__
60 #define __FUNCT__ "PCBDDCReuseMumpsReset"
61 static PetscErrorCode PCBDDCReuseMumpsReset(PCBDDCReuseMumps reuse)
62 {
63   PetscErrorCode ierr;
64 
65   PetscFunctionBegin;
66   ierr = MatDestroy(&reuse->F);CHKERRQ(ierr);
67   ierr = VecDestroy(&reuse->sol);CHKERRQ(ierr);
68   ierr = VecDestroy(&reuse->rhs);CHKERRQ(ierr);
69   ierr = PCDestroy(&reuse->interior_solver);CHKERRQ(ierr);
70   ierr = PCDestroy(&reuse->correction_solver);CHKERRQ(ierr);
71   ierr = ISDestroy(&reuse->is_R);CHKERRQ(ierr);
72   ierr = ISDestroy(&reuse->is_B);CHKERRQ(ierr);
73   ierr = VecScatterDestroy(&reuse->correction_scatter_B);CHKERRQ(ierr);
74   ierr = VecDestroy(&reuse->sol_B);CHKERRQ(ierr);
75   ierr = VecDestroy(&reuse->rhs_B);CHKERRQ(ierr);
76   PetscFunctionReturn(0);
77 }
78 
79 #undef __FUNCT__
80 #define __FUNCT__ "PCBDDCMumpsInteriorSolve_Private"
81 static PetscErrorCode PCBDDCMumpsInteriorSolve_Private(PC pc, Vec rhs, Vec sol, PetscBool transpose)
82 {
83   PCBDDCReuseMumps ctx;
84   PetscScalar      *array,*array_mumps;
85 #if defined(PETSC_HAVE_MUMPS)
86   PetscInt         ival;
87 #endif
88   PetscErrorCode   ierr;
89 
90   PetscFunctionBegin;
91   ierr = PCShellGetContext(pc,(void **)&ctx);CHKERRQ(ierr);
92 #if defined(PETSC_HAVE_MUMPS)
93   ierr = MatMumpsGetIcntl(ctx->F,26,&ival);CHKERRQ(ierr);
94   ierr = MatMumpsSetIcntl(ctx->F,26,0);CHKERRQ(ierr);
95 #endif
96   /* copy rhs into factored matrix workspace (can it be avoided?, MatSolve_MUMPS has another copy b->x internally) */
97   ierr = VecGetArrayRead(rhs,(const PetscScalar**)&array);CHKERRQ(ierr);
98   ierr = VecGetArray(ctx->rhs,&array_mumps);CHKERRQ(ierr);
99   ierr = PetscMemcpy(array_mumps,array,ctx->n_I*sizeof(PetscScalar));CHKERRQ(ierr);
100   ierr = VecRestoreArray(ctx->rhs,&array_mumps);CHKERRQ(ierr);
101   ierr = VecRestoreArrayRead(rhs,(const PetscScalar**)&array);CHKERRQ(ierr);
102 
103   if (transpose) {
104     ierr = MatSolveTranspose(ctx->F,ctx->rhs,ctx->sol);CHKERRQ(ierr);
105   } else {
106     ierr = MatSolve(ctx->F,ctx->rhs,ctx->sol);CHKERRQ(ierr);
107   }
108 
109   /* get back data to caller worskpace */
110   ierr = VecGetArrayRead(ctx->sol,(const PetscScalar**)&array_mumps);CHKERRQ(ierr);
111   ierr = VecGetArray(sol,&array);CHKERRQ(ierr);
112   ierr = PetscMemcpy(array,array_mumps,ctx->n_I*sizeof(PetscScalar));CHKERRQ(ierr);
113   ierr = VecRestoreArray(sol,&array);CHKERRQ(ierr);
114   ierr = VecRestoreArrayRead(ctx->sol,(const PetscScalar**)&array_mumps);CHKERRQ(ierr);
115 #if defined(PETSC_HAVE_MUMPS)
116   ierr = MatMumpsSetIcntl(ctx->F,26,ival);CHKERRQ(ierr);
117 #endif
118   PetscFunctionReturn(0);
119 }
120 
121 #undef __FUNCT__
122 #define __FUNCT__ "PCBDDCMumpsInteriorSolve"
123 static PetscErrorCode PCBDDCMumpsInteriorSolve(PC pc, Vec rhs, Vec sol)
124 {
125   PetscErrorCode   ierr;
126 
127   PetscFunctionBegin;
128   ierr = PCBDDCMumpsInteriorSolve_Private(pc,rhs,sol,PETSC_FALSE);CHKERRQ(ierr);
129   PetscFunctionReturn(0);
130 }
131 
132 #undef __FUNCT__
133 #define __FUNCT__ "PCBDDCMumpsInteriorSolveTranspose"
134 static PetscErrorCode PCBDDCMumpsInteriorSolveTranspose(PC pc, Vec rhs, Vec sol)
135 {
136   PetscErrorCode   ierr;
137 
138   PetscFunctionBegin;
139   ierr = PCBDDCMumpsInteriorSolve_Private(pc,rhs,sol,PETSC_TRUE);CHKERRQ(ierr);
140   PetscFunctionReturn(0);
141 }
142 
143 #undef __FUNCT__
144 #define __FUNCT__ "PCBDDCComputeExplicitSchur"
145 static PetscErrorCode PCBDDCComputeExplicitSchur(Mat M, PetscBool issym, MatReuse reuse, Mat *S)
146 {
147   Mat            B, C, D, Bd, Cd, AinvBd;
148   KSP            ksp;
149   PC             pc;
150   PetscBool      isLU, isILU, isCHOL, Bdense, Cdense;
151   PetscReal      fill = 2.0;
152   PetscInt       n_I;
153   PetscMPIInt    size;
154   PetscErrorCode ierr;
155 
156   PetscFunctionBegin;
157   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)M),&size);CHKERRQ(ierr);
158   if (size != 1) {
159     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not for parallel matrices");
160   }
161   if (reuse == MAT_REUSE_MATRIX) {
162     PetscBool Sdense;
163 
164     ierr = PetscObjectTypeCompare((PetscObject)*S, MATSEQDENSE, &Sdense);CHKERRQ(ierr);
165     if (!Sdense) {
166       SETERRQ(PetscObjectComm((PetscObject)M),PETSC_ERR_SUP,"S should dense");
167     }
168   }
169   ierr = MatSchurComplementGetSubMatrices(M, NULL, NULL, &B, &C, &D);CHKERRQ(ierr);
170   ierr = MatSchurComplementGetKSP(M, &ksp);CHKERRQ(ierr);
171   ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr);
172   ierr = PetscObjectTypeCompare((PetscObject) pc, PCLU, &isLU);CHKERRQ(ierr);
173   ierr = PetscObjectTypeCompare((PetscObject) pc, PCILU, &isILU);CHKERRQ(ierr);
174   ierr = PetscObjectTypeCompare((PetscObject) pc, PCCHOLESKY, &isCHOL);CHKERRQ(ierr);
175   ierr = PetscObjectTypeCompare((PetscObject) B, MATSEQDENSE, &Bdense);CHKERRQ(ierr);
176   ierr = PetscObjectTypeCompare((PetscObject) C, MATSEQDENSE, &Cdense);CHKERRQ(ierr);
177   ierr = MatGetSize(B,&n_I,NULL);CHKERRQ(ierr);
178   if (n_I) {
179     if (!Bdense) {
180       ierr = MatConvert(B, MATSEQDENSE, MAT_INITIAL_MATRIX, &Bd);CHKERRQ(ierr);
181     } else {
182       Bd = B;
183     }
184 
185     if (isLU || isILU || isCHOL) {
186       Mat fact;
187       ierr = KSPSetUp(ksp);CHKERRQ(ierr);
188       ierr = PCFactorGetMatrix(pc, &fact);CHKERRQ(ierr);
189       ierr = MatDuplicate(Bd, MAT_DO_NOT_COPY_VALUES, &AinvBd);CHKERRQ(ierr);
190       ierr = MatMatSolve(fact, Bd, AinvBd);CHKERRQ(ierr);
191     } else {
192       PetscBool ex = PETSC_TRUE;
193 
194       if (ex) {
195         Mat Ainvd;
196 
197         ierr = PCComputeExplicitOperator(pc, &Ainvd);CHKERRQ(ierr);
198         ierr = MatMatMult(Ainvd, Bd, MAT_INITIAL_MATRIX, fill, &AinvBd);CHKERRQ(ierr);
199         ierr = MatDestroy(&Ainvd);CHKERRQ(ierr);
200       } else {
201         Vec         sol,rhs;
202         PetscScalar *arrayrhs,*arraysol;
203         PetscInt    i,nrhs,n;
204 
205         ierr = MatDuplicate(Bd, MAT_DO_NOT_COPY_VALUES, &AinvBd);CHKERRQ(ierr);
206         ierr = MatGetSize(Bd,&n,&nrhs);CHKERRQ(ierr);
207         ierr = MatDenseGetArray(Bd,&arrayrhs);CHKERRQ(ierr);
208         ierr = MatDenseGetArray(AinvBd,&arraysol);CHKERRQ(ierr);
209         ierr = KSPGetSolution(ksp,&sol);CHKERRQ(ierr);
210         ierr = KSPGetRhs(ksp,&rhs);CHKERRQ(ierr);
211         for (i=0;i<nrhs;i++) {
212           ierr = VecPlaceArray(rhs,arrayrhs+i*n);CHKERRQ(ierr);
213           ierr = VecPlaceArray(sol,arraysol+i*n);CHKERRQ(ierr);
214           ierr = KSPSolve(ksp,rhs,sol);CHKERRQ(ierr);
215           ierr = VecResetArray(rhs);CHKERRQ(ierr);
216           ierr = VecResetArray(sol);CHKERRQ(ierr);
217         }
218         ierr = MatDenseRestoreArray(Bd,&arrayrhs);CHKERRQ(ierr);
219         ierr = MatDenseRestoreArray(AinvBd,&arrayrhs);CHKERRQ(ierr);
220       }
221     }
222     if (!Bdense & !issym) {
223       ierr = MatDestroy(&Bd);CHKERRQ(ierr);
224     }
225 
226     if (!issym) {
227       if (!Cdense) {
228         ierr = MatConvert(C, MATSEQDENSE, MAT_INITIAL_MATRIX, &Cd);CHKERRQ(ierr);
229       } else {
230         Cd = C;
231       }
232       ierr = MatMatMult(Cd, AinvBd, reuse, fill, S);CHKERRQ(ierr);
233       if (!Cdense) {
234         ierr = MatDestroy(&Cd);CHKERRQ(ierr);
235       }
236     } else {
237       ierr = MatTransposeMatMult(Bd, AinvBd, reuse, fill, S);CHKERRQ(ierr);
238       if (!Bdense) {
239         ierr = MatDestroy(&Bd);CHKERRQ(ierr);
240       }
241     }
242     ierr = MatDestroy(&AinvBd);CHKERRQ(ierr);
243   }
244 
245   if (D) {
246     Mat       Dd;
247     PetscBool Ddense;
248 
249     ierr = PetscObjectTypeCompare((PetscObject)D,MATSEQDENSE,&Ddense);CHKERRQ(ierr);
250     if (!Ddense) {
251       ierr = MatConvert(D, MATSEQDENSE, MAT_INITIAL_MATRIX, &Dd);CHKERRQ(ierr);
252     } else {
253       Dd = D;
254     }
255     if (n_I) {
256       ierr = MatAYPX(*S,-1.0,Dd,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
257     } else {
258       if (reuse == MAT_INITIAL_MATRIX) {
259         ierr = MatDuplicate(Dd,MAT_COPY_VALUES,S);CHKERRQ(ierr);
260       } else {
261         ierr = MatCopy(Dd,*S,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
262       }
263     }
264     if (!Ddense) {
265       ierr = MatDestroy(&Dd);CHKERRQ(ierr);
266     }
267   } else {
268     ierr = MatScale(*S,-1.0);CHKERRQ(ierr);
269   }
270   PetscFunctionReturn(0);
271 }
272 
273 #undef __FUNCT__
274 #define __FUNCT__ "PCBDDCSubSchursSetUp"
275 PetscErrorCode PCBDDCSubSchursSetUp(PCBDDCSubSchurs sub_schurs, Mat Ain, Mat Sin, PetscInt xadj[], PetscInt adjncy[], PetscInt nlayers, PetscBool faster_deluxe, PetscBool compute_Stilda, PetscBool reuse_solvers,PetscBool use_edges, PetscBool use_faces)
276 {
277   Mat                    F,A_II,A_IB,A_BI,A_BB,AE_II;
278   Mat                    S_all;
279   Mat                    global_schur_subsets,work_mat;
280   Mat                    S_Ej_tilda_all,S_Ej_inv_all;
281   ISLocalToGlobalMapping l2gmap_subsets;
282   IS                     is_I,is_I_layer,temp_is;
283   PetscInt               *nnz,*all_local_idx_G,*all_local_idx_N;
284   PetscInt               *auxnum1,*auxnum2,*all_local_idx_G_rep;
285   PetscInt               i,subset_size,max_subset_size;
286   PetscInt               extra,local_size,global_size;
287   PetscBLASInt           B_N,B_ierr,B_lwork,*pivots;
288   PetscScalar            *Bwork;
289   PetscSubcomm           subcomm;
290   PetscMPIInt            color,rank;
291   MPI_Comm               comm_n;
292   PetscErrorCode         ierr;
293 
294   PetscFunctionBegin;
295   /* update info in sub_schurs */
296   ierr = MatDestroy(&sub_schurs->A);CHKERRQ(ierr);
297   ierr = MatDestroy(&sub_schurs->S);CHKERRQ(ierr);
298   if (Ain) {
299     PetscBool isseqaij;
300 
301     ierr = PetscObjectTypeCompare((PetscObject)Ain,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
302     if (isseqaij) {
303       ierr = PetscObjectReference((PetscObject)Ain);CHKERRQ(ierr);
304       sub_schurs->A = Ain;
305     } else { /* SeqBAIJ matrices does not support symmetry checking, SeqSBAIJ does not support MatPermute */
306       ierr = MatConvert(Ain,MATSEQAIJ,MAT_INITIAL_MATRIX,&sub_schurs->A);CHKERRQ(ierr);
307     }
308   }
309   ierr = PetscObjectReference((PetscObject)Sin);CHKERRQ(ierr);
310   sub_schurs->S = Sin;
311   if (sub_schurs->use_mumps) {
312     sub_schurs->use_mumps = (PetscBool)(!!sub_schurs->A);
313   }
314 
315   /* preliminary checks */
316   if (!sub_schurs->use_mumps && compute_Stilda) {
317     SETERRQ(PetscObjectComm((PetscObject)sub_schurs->l2gmap),PETSC_ERR_SUP,"Adaptive selection of constraints requires MUMPS");
318   }
319   /* determine if we are dealing with hermitian positive definite problems */
320   sub_schurs->is_hermitian = PETSC_FALSE;
321   sub_schurs->is_posdef = PETSC_FALSE;
322   if (sub_schurs->A) {
323     PetscInt lsize;
324 
325     ierr = MatGetSize(sub_schurs->A,&lsize,NULL);CHKERRQ(ierr);
326     if (lsize) {
327       ierr = MatIsHermitian(sub_schurs->A,0.0,&sub_schurs->is_hermitian);CHKERRQ(ierr);
328       if (sub_schurs->is_hermitian) {
329         PetscScalar val;
330         Vec         vec1,vec2;
331 
332         ierr = MatCreateVecs(sub_schurs->A,&vec1,&vec2);CHKERRQ(ierr);
333         ierr = VecSetRandom(vec1,NULL);
334         ierr = VecCopy(vec1,vec2);CHKERRQ(ierr);
335         ierr = MatMult(sub_schurs->A,vec2,vec1);CHKERRQ(ierr);
336         ierr = VecDot(vec1,vec2,&val);CHKERRQ(ierr);
337         if (PetscRealPart(val) > 0. && PetscAbsReal(PetscImaginaryPart(val)) < PETSC_SMALL) sub_schurs->is_posdef = PETSC_TRUE;
338         ierr = VecDestroy(&vec1);CHKERRQ(ierr);
339         ierr = VecDestroy(&vec2);CHKERRQ(ierr);
340       }
341     } else {
342       sub_schurs->is_hermitian = PETSC_TRUE;
343       sub_schurs->is_posdef = PETSC_TRUE;
344     }
345     if (compute_Stilda && (!sub_schurs->is_hermitian || !sub_schurs->is_posdef)) {
346       SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"General matrix pencils are not currently supported (%d,%d)",sub_schurs->is_hermitian,sub_schurs->is_posdef);
347     }
348   }
349   /* restrict work on active processes */
350   color = 0;
351   if (!sub_schurs->n_subs) color = 1; /* this can happen if we are in a multilevel case or if the subdomain is disconnected */
352   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&rank);CHKERRQ(ierr);
353   ierr = PetscSubcommCreate(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&subcomm);CHKERRQ(ierr);
354   ierr = PetscSubcommSetNumber(subcomm,2);CHKERRQ(ierr);
355   ierr = PetscSubcommSetTypeGeneral(subcomm,color,rank);CHKERRQ(ierr);
356   ierr = PetscCommDuplicate(PetscSubcommChild(subcomm),&comm_n,NULL);CHKERRQ(ierr);
357   ierr = PetscSubcommDestroy(&subcomm);CHKERRQ(ierr);
358   if (!sub_schurs->n_subs) {
359     ierr = PetscCommDestroy(&comm_n);CHKERRQ(ierr);
360     PetscFunctionReturn(0);
361   }
362 
363   /* get Schur complement matrices */
364   if (!sub_schurs->use_mumps) {
365     Mat       tA_IB,tA_BI,tA_BB;
366     PetscBool isseqaij;
367     ierr = MatSchurComplementGetSubMatrices(sub_schurs->S,&A_II,NULL,&tA_IB,&tA_BI,&tA_BB);CHKERRQ(ierr);
368     ierr = PetscObjectTypeCompare((PetscObject)tA_BB,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
369     if (!isseqaij) {
370       ierr = MatConvert(tA_BB,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_BB);CHKERRQ(ierr);
371       ierr = MatConvert(tA_IB,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_IB);CHKERRQ(ierr);
372       ierr = MatConvert(tA_BI,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_BI);CHKERRQ(ierr);
373     } else {
374       ierr = PetscObjectReference((PetscObject)tA_BB);CHKERRQ(ierr);
375       A_BB = tA_BB;
376       ierr = PetscObjectReference((PetscObject)tA_IB);CHKERRQ(ierr);
377       A_IB = tA_IB;
378       ierr = PetscObjectReference((PetscObject)tA_BI);CHKERRQ(ierr);
379       A_BI = tA_BI;
380     }
381   } else {
382     A_II = NULL;
383     A_IB = NULL;
384     A_BI = NULL;
385     A_BB = NULL;
386   }
387   S_all = NULL;
388   S_Ej_tilda_all = NULL;
389   S_Ej_inv_all = NULL;
390 
391   /* determine interior problems */
392   ierr = ISGetLocalSize(sub_schurs->is_I,&i);CHKERRQ(ierr);
393   if (nlayers >= 0 && i) { /* Interior problems can be different from the original one */
394     PetscBT                touched;
395     const PetscInt*        idx_B;
396     PetscInt               n_I,n_B,n_local_dofs,n_prev_added,j,layer,*local_numbering;
397 
398     if (xadj == NULL || adjncy == NULL) {
399       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cannot request layering without adjacency");
400     }
401     /* get sizes */
402     ierr = ISGetLocalSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr);
403     ierr = ISGetLocalSize(sub_schurs->is_B,&n_B);CHKERRQ(ierr);
404 
405     ierr = PetscMalloc1(n_I+n_B,&local_numbering);CHKERRQ(ierr);
406     ierr = PetscBTCreate(n_I+n_B,&touched);CHKERRQ(ierr);
407     ierr = PetscBTMemzero(n_I+n_B,touched);CHKERRQ(ierr);
408 
409     /* all boundary dofs must be skipped when adding layers */
410     ierr = ISGetIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr);
411     for (j=0;j<n_B;j++) {
412       ierr = PetscBTSet(touched,idx_B[j]);CHKERRQ(ierr);
413     }
414     ierr = PetscMemcpy(local_numbering,idx_B,n_B*sizeof(PetscInt));CHKERRQ(ierr);
415     ierr = ISRestoreIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr);
416 
417     /* add prescribed number of layers of dofs */
418     n_local_dofs = n_B;
419     n_prev_added = n_B;
420     for (layer=0;layer<nlayers;layer++) {
421       PetscInt n_added;
422       if (n_local_dofs == n_I+n_B) break;
423       if (n_local_dofs > n_I+n_B) {
424         SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error querying layer %d. Out of bound access (%d > %d)",layer,n_local_dofs,n_I+n_B);
425       }
426       ierr = PCBDDCAdjGetNextLayer_Private(local_numbering+n_local_dofs,n_prev_added,touched,xadj,adjncy,&n_added);CHKERRQ(ierr);
427       n_prev_added = n_added;
428       n_local_dofs += n_added;
429       if (!n_added) break;
430     }
431     ierr = PetscBTDestroy(&touched);CHKERRQ(ierr);
432 
433     /* IS for I layer dofs in original numbering */
434     ierr = ISCreateGeneral(PetscObjectComm((PetscObject)sub_schurs->is_I),n_local_dofs-n_B,local_numbering+n_B,PETSC_COPY_VALUES,&is_I_layer);CHKERRQ(ierr);
435     ierr = PetscFree(local_numbering);CHKERRQ(ierr);
436     ierr = ISSort(is_I_layer);CHKERRQ(ierr);
437     /* IS for I layer dofs in I numbering */
438     if (!sub_schurs->use_mumps) {
439       ISLocalToGlobalMapping ItoNmap;
440       ierr = ISLocalToGlobalMappingCreateIS(sub_schurs->is_I,&ItoNmap);CHKERRQ(ierr);
441       ierr = ISGlobalToLocalMappingApplyIS(ItoNmap,IS_GTOLM_DROP,is_I_layer,&is_I);CHKERRQ(ierr);
442       ierr = ISLocalToGlobalMappingDestroy(&ItoNmap);CHKERRQ(ierr);
443 
444       /* II block */
445       ierr = MatGetSubMatrix(A_II,is_I,is_I,MAT_INITIAL_MATRIX,&AE_II);CHKERRQ(ierr);
446     }
447   } else {
448     PetscInt n_I;
449 
450     /* IS for I dofs in original numbering */
451     ierr = PetscObjectReference((PetscObject)sub_schurs->is_I);CHKERRQ(ierr);
452     is_I_layer = sub_schurs->is_I;
453 
454     /* IS for I dofs in I numbering (strided 1) */
455     if (!sub_schurs->use_mumps) {
456       ierr = ISGetSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr);
457       ierr = ISCreateStride(PetscObjectComm((PetscObject)sub_schurs->is_I),n_I,0,1,&is_I);CHKERRQ(ierr);
458 
459       /* II block is the same */
460       ierr = PetscObjectReference((PetscObject)A_II);CHKERRQ(ierr);
461       AE_II = A_II;
462     }
463   }
464 
465   /* Get info on subset sizes and sum of all subsets sizes */
466   max_subset_size = 0;
467   local_size = 0;
468   for (i=0;i<sub_schurs->n_subs;i++) {
469     ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
470     max_subset_size = PetscMax(subset_size,max_subset_size);
471     local_size += subset_size;
472   }
473 
474   /* Work arrays for local indices */
475   extra = 0;
476   if (sub_schurs->use_mumps) {
477     ierr = ISGetLocalSize(is_I_layer,&extra);CHKERRQ(ierr);
478   }
479   ierr = PetscMalloc1(local_size+extra,&all_local_idx_N);CHKERRQ(ierr);
480   if (extra) {
481     const PetscInt *idxs;
482     ierr = ISGetIndices(is_I_layer,&idxs);CHKERRQ(ierr);
483     ierr = PetscMemcpy(all_local_idx_N,idxs,extra*sizeof(PetscInt));CHKERRQ(ierr);
484     ierr = ISRestoreIndices(is_I_layer,&idxs);CHKERRQ(ierr);
485   }
486   ierr = PetscMalloc1(local_size,&nnz);CHKERRQ(ierr);
487   ierr = PetscMalloc2(sub_schurs->n_subs,&auxnum1,sub_schurs->n_subs,&auxnum2);CHKERRQ(ierr);
488 
489   /* Get local indices in local numbering */
490   local_size = 0;
491   for (i=0;i<sub_schurs->n_subs;i++) {
492     PetscInt j;
493     const    PetscInt *idxs;
494 
495     ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
496     ierr = ISGetIndices(sub_schurs->is_subs[i],&idxs);CHKERRQ(ierr);
497     /* start (smallest in global ordering) and multiplicity */
498     auxnum1[i] = idxs[0];
499     auxnum2[i] = subset_size;
500     /* subset indices in local numbering */
501     ierr = PetscMemcpy(all_local_idx_N+local_size+extra,idxs,subset_size*sizeof(PetscInt));CHKERRQ(ierr);
502     ierr = ISRestoreIndices(sub_schurs->is_subs[i],&idxs);CHKERRQ(ierr);
503     for (j=0;j<subset_size;j++) nnz[local_size+j] = subset_size;
504     local_size += subset_size;
505   }
506 
507   /* allocate extra workspace needed only for GETRI */
508   Bwork = NULL;
509   pivots = NULL;
510   if (sub_schurs->n_subs && !sub_schurs->is_hermitian) {
511     PetscScalar lwork;
512 
513     B_lwork = -1;
514     ierr = PetscBLASIntCast(local_size,&B_N);CHKERRQ(ierr);
515     ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
516     PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,Bwork,&B_N,pivots,&lwork,&B_lwork,&B_ierr));
517     ierr = PetscFPTrapPop();CHKERRQ(ierr);
518     if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to GETRI Lapack routine %d",(int)B_ierr);
519     ierr = PetscBLASIntCast((PetscInt)PetscRealPart(lwork),&B_lwork);CHKERRQ(ierr);
520     ierr = PetscMalloc2(B_lwork,&Bwork,B_N,&pivots);CHKERRQ(ierr);
521   }
522 
523   /* prepare parallel matrices for summing up properly schurs on subsets */
524   ierr = PCBDDCSubsetNumbering(comm_n,sub_schurs->l2gmap,sub_schurs->n_subs,auxnum1,auxnum2,&global_size,&all_local_idx_G_rep);CHKERRQ(ierr);
525   ierr = PetscMalloc1(local_size,&all_local_idx_G);CHKERRQ(ierr);
526   local_size = 0;
527   for (i=0;i<sub_schurs->n_subs;i++) {
528     PetscInt j;
529     for (j=0;j<auxnum2[i];j++) all_local_idx_G[local_size++] = all_local_idx_G_rep[i] + j;
530   }
531   ierr = PetscFree(all_local_idx_G_rep);CHKERRQ(ierr);
532   ierr = PetscFree2(auxnum1,auxnum2);CHKERRQ(ierr);
533   ierr = ISLocalToGlobalMappingCreate(comm_n,1,local_size,all_local_idx_G,PETSC_COPY_VALUES,&l2gmap_subsets);CHKERRQ(ierr);
534   ierr = MatCreateIS(comm_n,1,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size,l2gmap_subsets,&work_mat);CHKERRQ(ierr);
535   ierr = ISLocalToGlobalMappingDestroy(&l2gmap_subsets);CHKERRQ(ierr);
536   ierr = MatCreate(PetscObjectComm((PetscObject)work_mat),&global_schur_subsets);CHKERRQ(ierr);
537   ierr = MatSetSizes(global_schur_subsets,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size);CHKERRQ(ierr);
538   ierr = MatSetType(global_schur_subsets,MATMPIAIJ);CHKERRQ(ierr);
539 
540   /* subset indices in local boundary numbering */
541   if (!sub_schurs->is_Ej_all) {
542     PetscInt *all_local_idx_B;
543 
544     ierr = PetscMalloc1(local_size,&all_local_idx_B);CHKERRQ(ierr);
545     ierr = ISGlobalToLocalMappingApply(sub_schurs->BtoNmap,IS_GTOLM_DROP,local_size,all_local_idx_N+extra,&subset_size,all_local_idx_B);CHKERRQ(ierr);
546     if (subset_size != local_size) {
547       SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in sub_schurs serial (BtoNmap)! %d != %d\n",subset_size,local_size);
548     }
549     ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size,all_local_idx_B,PETSC_OWN_POINTER,&sub_schurs->is_Ej_all);CHKERRQ(ierr);
550   }
551 
552   /* Local matrix of all local Schur on subsets (transposed) */
553   if (!sub_schurs->S_Ej_all) {
554     ierr = MatCreate(PETSC_COMM_SELF,&sub_schurs->S_Ej_all);CHKERRQ(ierr);
555     ierr = MatSetSizes(sub_schurs->S_Ej_all,PETSC_DECIDE,PETSC_DECIDE,local_size,local_size);CHKERRQ(ierr);
556     ierr = MatSetType(sub_schurs->S_Ej_all,MATAIJ);CHKERRQ(ierr);
557     ierr = MatSeqAIJSetPreallocation(sub_schurs->S_Ej_all,0,nnz);CHKERRQ(ierr);
558   } else {
559     ierr = MatZeroEntries(sub_schurs->S_Ej_all);CHKERRQ(ierr);
560   }
561 
562   /* Compute Schur complements explicitly */
563   F = NULL;
564   ierr = PetscBTMemzero(sub_schurs->n_subs,sub_schurs->computed_Stilda_subs);CHKERRQ(ierr);
565   if (!sub_schurs->use_mumps) {
566     Mat         S_Ej_expl;
567     PetscScalar *work;
568     PetscInt    j,*dummy_idx;
569     PetscBool   Sdense;
570 
571     ierr = PetscMalloc2(max_subset_size,&dummy_idx,max_subset_size*max_subset_size,&work);CHKERRQ(ierr);
572     local_size = 0;
573     for (i=0;i<sub_schurs->n_subs;i++) {
574       IS  is_subset_B;
575       Mat AE_EE,AE_IE,AE_EI,S_Ej;
576 
577       /* subsets in original and boundary numbering */
578       ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,sub_schurs->is_subs[i],&is_subset_B);CHKERRQ(ierr);
579       /* EE block */
580       ierr = MatGetSubMatrix(A_BB,is_subset_B,is_subset_B,MAT_INITIAL_MATRIX,&AE_EE);CHKERRQ(ierr);
581       /* IE block */
582       ierr = MatGetSubMatrix(A_IB,is_I,is_subset_B,MAT_INITIAL_MATRIX,&AE_IE);CHKERRQ(ierr);
583       /* EI block */
584       if (sub_schurs->is_hermitian) {
585         ierr = MatCreateTranspose(AE_IE,&AE_EI);CHKERRQ(ierr);
586       } else {
587         ierr = MatGetSubMatrix(A_BI,is_subset_B,is_I,MAT_INITIAL_MATRIX,&AE_EI);CHKERRQ(ierr);
588       }
589       ierr = ISDestroy(&is_subset_B);CHKERRQ(ierr);
590       ierr = MatCreateSchurComplement(AE_II,AE_II,AE_IE,AE_EI,AE_EE,&S_Ej);CHKERRQ(ierr);
591       ierr = MatDestroy(&AE_EE);CHKERRQ(ierr);
592       ierr = MatDestroy(&AE_IE);CHKERRQ(ierr);
593       ierr = MatDestroy(&AE_EI);CHKERRQ(ierr);
594       if (AE_II == A_II) { /* we can reuse the same ksp */
595         KSP ksp;
596         ierr = MatSchurComplementGetKSP(sub_schurs->S,&ksp);CHKERRQ(ierr);
597         ierr = MatSchurComplementSetKSP(S_Ej,ksp);CHKERRQ(ierr);
598       } else { /* build new ksp object which inherits ksp and pc types from the original one */
599         KSP       origksp,schurksp;
600         PC        origpc,schurpc;
601         KSPType   ksp_type;
602         PetscInt  n_internal;
603         PetscBool ispcnone;
604 
605         ierr = MatSchurComplementGetKSP(sub_schurs->S,&origksp);CHKERRQ(ierr);
606         ierr = MatSchurComplementGetKSP(S_Ej,&schurksp);CHKERRQ(ierr);
607         ierr = KSPGetType(origksp,&ksp_type);CHKERRQ(ierr);
608         ierr = KSPSetType(schurksp,ksp_type);CHKERRQ(ierr);
609         ierr = KSPGetPC(schurksp,&schurpc);CHKERRQ(ierr);
610         ierr = KSPGetPC(origksp,&origpc);CHKERRQ(ierr);
611         ierr = PetscObjectTypeCompare((PetscObject)origpc,PCNONE,&ispcnone);CHKERRQ(ierr);
612         if (!ispcnone) {
613           PCType pc_type;
614           ierr = PCGetType(origpc,&pc_type);CHKERRQ(ierr);
615           ierr = PCSetType(schurpc,pc_type);CHKERRQ(ierr);
616         } else {
617           ierr = PCSetType(schurpc,PCLU);CHKERRQ(ierr);
618         }
619         ierr = ISGetSize(is_I,&n_internal);CHKERRQ(ierr);
620         if (n_internal) { /* UMFPACK gives error with 0 sized problems */
621           MatSolverPackage solver=NULL;
622           ierr = PCFactorGetMatSolverPackage(origpc,(const MatSolverPackage*)&solver);CHKERRQ(ierr);
623           if (solver) {
624             ierr = PCFactorSetMatSolverPackage(schurpc,solver);CHKERRQ(ierr);
625           }
626         }
627         ierr = KSPSetUp(schurksp);CHKERRQ(ierr);
628       }
629       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
630       ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&S_Ej_expl);CHKERRQ(ierr);
631       ierr = PCBDDCComputeExplicitSchur(S_Ej,sub_schurs->is_hermitian,MAT_REUSE_MATRIX,&S_Ej_expl);CHKERRQ(ierr);
632       ierr = PetscObjectTypeCompare((PetscObject)S_Ej_expl,MATSEQDENSE,&Sdense);CHKERRQ(ierr);
633       if (Sdense) {
634         for (j=0;j<subset_size;j++) {
635           dummy_idx[j]=local_size+j;
636         }
637         ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
638       } else {
639         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not yet implemented for sparse matrices");
640       }
641       ierr = MatDestroy(&S_Ej);CHKERRQ(ierr);
642       ierr = MatDestroy(&S_Ej_expl);CHKERRQ(ierr);
643       local_size += subset_size;
644     }
645     ierr = PetscFree2(dummy_idx,work);CHKERRQ(ierr);
646     /* free */
647     ierr = ISDestroy(&is_I);CHKERRQ(ierr);
648     ierr = MatDestroy(&AE_II);CHKERRQ(ierr);
649     ierr = PetscFree(all_local_idx_N);CHKERRQ(ierr);
650   } else {
651     Mat         A;
652     IS          is_A_all;
653     PetscScalar *work,*S_data;
654     PetscInt    *idxs_schur,n_I,n_I_all,*dummy_idx,size_schur,size_active_schur,cum,cum2;
655     PetscBool   mumps_S;
656 
657     /* get working mat */
658     ierr = ISGetLocalSize(is_I_layer,&n_I);CHKERRQ(ierr);
659     if (!sub_schurs->is_dir) {
660       ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size+n_I,all_local_idx_N,PETSC_COPY_VALUES,&is_A_all);CHKERRQ(ierr);
661       size_schur = local_size;
662     } else {
663       IS list[2];
664 
665       ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size+n_I,all_local_idx_N,PETSC_COPY_VALUES,&list[0]);CHKERRQ(ierr);
666       list[1] = sub_schurs->is_dir;
667       ierr = ISConcatenate(PETSC_COMM_SELF,2,list,&is_A_all);CHKERRQ(ierr);
668       ierr = ISDestroy(&list[0]);CHKERRQ(ierr);
669       ierr = ISGetLocalSize(sub_schurs->is_dir,&size_schur);CHKERRQ(ierr);
670       size_schur += local_size;
671     }
672     ierr = PetscFree(all_local_idx_N);CHKERRQ(ierr);
673     size_active_schur = local_size; /* size active schurs does not count any dirichlet dof on the interface */
674     ierr = MatGetSubMatrix(sub_schurs->A,is_A_all,is_A_all,MAT_INITIAL_MATRIX,&A);CHKERRQ(ierr);
675     ierr = MatSetOptionsPrefix(A,"sub_schurs_");CHKERRQ(ierr);
676 
677     if (n_I) {
678       if (sub_schurs->is_hermitian && sub_schurs->is_posdef) {
679         ierr = MatSetOption(A,MAT_SPD,PETSC_TRUE);CHKERRQ(ierr);
680         ierr = MatGetFactor(A,MATSOLVERMUMPS,MAT_FACTOR_CHOLESKY,&F);CHKERRQ(ierr);
681       } else {
682         ierr = MatSetOption(A,MAT_SYMMETRIC,sub_schurs->is_hermitian);CHKERRQ(ierr);
683         ierr = MatGetFactor(A,MATSOLVERMUMPS,MAT_FACTOR_LU,&F);CHKERRQ(ierr);
684       }
685 
686       /* subsets ordered last */
687       ierr = PetscMalloc1(size_schur,&idxs_schur);CHKERRQ(ierr);
688       for (i=0;i<size_schur;i++) {
689         idxs_schur[i] = n_I+i+1;
690       }
691 #if defined(PETSC_HAVE_MUMPS)
692       ierr = MatMumpsSetSchurIndices(F,size_schur,idxs_schur);CHKERRQ(ierr);
693 #endif
694       ierr = PetscFree(idxs_schur);CHKERRQ(ierr);
695 
696       /* factorization step */
697       if (sub_schurs->is_hermitian && sub_schurs->is_posdef) {
698         ierr = MatCholeskyFactorSymbolic(F,A,NULL,NULL);CHKERRQ(ierr);
699 #if defined(PETSC_HAVE_MUMPS) /* be sure that icntl 19 is not set by command line */
700         ierr = MatMumpsSetIcntl(F,19,2);CHKERRQ(ierr);
701 #endif
702         ierr = MatCholeskyFactorNumeric(F,A,NULL);CHKERRQ(ierr);
703       } else {
704         ierr = MatLUFactorSymbolic(F,A,NULL,NULL,NULL);CHKERRQ(ierr);
705 #if defined(PETSC_HAVE_MUMPS) /* be sure that icntl 19 is not set by command line */
706         ierr = MatMumpsSetIcntl(F,19,3);CHKERRQ(ierr);
707 #endif
708         ierr = MatLUFactorNumeric(F,A,NULL);CHKERRQ(ierr);
709       }
710 
711       /* get explicit Schur Complement computed during numeric factorization */
712 #if defined(PETSC_HAVE_MUMPS)
713       ierr = MatMumpsGetSchurComplement(F,&S_all);CHKERRQ(ierr);
714 #endif
715 
716       /* we can reuse the solvers if we are not using the economic version */
717       ierr = ISGetLocalSize(sub_schurs->is_I,&n_I_all);CHKERRQ(ierr);
718       reuse_solvers = (PetscBool)(reuse_solvers && (n_I == n_I_all));
719       mumps_S = PETSC_TRUE;
720     } else { /* we can't use MUMPS when size_schur == size_of_the_problem */
721       ierr = MatConvert(A,MATSEQDENSE,MAT_INITIAL_MATRIX,&S_all);CHKERRQ(ierr);
722       reuse_solvers = PETSC_FALSE;
723       mumps_S = PETSC_FALSE;
724     }
725 
726     if (reuse_solvers) {
727       Mat              A_II;
728       Vec              vec1_B;
729       PCBDDCReuseMumps msolv_ctx;
730 
731       if (sub_schurs->reuse_mumps) {
732         ierr = PCBDDCReuseMumpsReset(sub_schurs->reuse_mumps);CHKERRQ(ierr);
733       } else {
734         ierr = PetscNew(&sub_schurs->reuse_mumps);CHKERRQ(ierr);
735       }
736       msolv_ctx = sub_schurs->reuse_mumps;
737       ierr = MatSchurComplementGetSubMatrices(sub_schurs->S,&A_II,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
738       ierr = MatGetSize(A_II,&msolv_ctx->n_I,NULL);CHKERRQ(ierr);
739       ierr = PetscObjectReference((PetscObject)F);CHKERRQ(ierr);
740       msolv_ctx->F = F;
741       ierr = MatCreateVecs(F,&msolv_ctx->sol,&msolv_ctx->rhs);CHKERRQ(ierr);
742 
743       /* interior solver */
744       ierr = PCCreate(PETSC_COMM_SELF,&msolv_ctx->interior_solver);CHKERRQ(ierr);
745       ierr = PCSetOperators(msolv_ctx->interior_solver,A_II,A_II);CHKERRQ(ierr);
746       ierr = PCSetType(msolv_ctx->interior_solver,PCSHELL);CHKERRQ(ierr);
747       ierr = PCShellSetContext(msolv_ctx->interior_solver,msolv_ctx);CHKERRQ(ierr);
748       ierr = PCShellSetApply(msolv_ctx->interior_solver,PCBDDCMumpsInteriorSolve);CHKERRQ(ierr);
749       ierr = PCShellSetApplyTranspose(msolv_ctx->interior_solver,PCBDDCMumpsInteriorSolveTranspose);CHKERRQ(ierr);
750 
751       /* correction solver */
752       ierr = PCCreate(PETSC_COMM_SELF,&msolv_ctx->correction_solver);CHKERRQ(ierr);
753       ierr = PCSetOperators(msolv_ctx->correction_solver,A,A);CHKERRQ(ierr);
754       ierr = PCSetType(msolv_ctx->correction_solver,PCSHELL);CHKERRQ(ierr);
755       ierr = PCShellSetContext(msolv_ctx->correction_solver,msolv_ctx);CHKERRQ(ierr);
756       ierr = PCShellSetApply(msolv_ctx->correction_solver,PCBDDCMumpsCorrectionSolve);CHKERRQ(ierr);
757       ierr = PCShellSetApplyTranspose(msolv_ctx->correction_solver,PCBDDCMumpsCorrectionSolveTranspose);CHKERRQ(ierr);
758 
759       /* scatter and vecs for Schur complement solver */
760       ierr = MatCreateVecs(S_all,&msolv_ctx->sol_B,&msolv_ctx->rhs_B);CHKERRQ(ierr);
761       ierr = MatCreateVecs(sub_schurs->S,&vec1_B,NULL);CHKERRQ(ierr);
762       ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,is_A_all,&msolv_ctx->is_B);CHKERRQ(ierr);
763       ierr = VecScatterCreate(vec1_B,msolv_ctx->is_B,msolv_ctx->sol_B,NULL,&msolv_ctx->correction_scatter_B);CHKERRQ(ierr);
764       ierr = VecDestroy(&vec1_B);CHKERRQ(ierr);
765       ierr = PetscObjectReference((PetscObject)is_A_all);CHKERRQ(ierr);
766       msolv_ctx->is_R = is_A_all;
767     }
768     ierr = MatDestroy(&A);CHKERRQ(ierr);
769     ierr = ISDestroy(&is_A_all);CHKERRQ(ierr);
770 
771     /* Work arrays */
772     ierr = PetscMalloc2(max_subset_size,&dummy_idx,max_subset_size*max_subset_size,&work);CHKERRQ(ierr);
773     if (compute_Stilda) {
774       ierr = MatCreate(PETSC_COMM_SELF,&S_Ej_tilda_all);CHKERRQ(ierr);
775       ierr = MatSetSizes(S_Ej_tilda_all,PETSC_DECIDE,PETSC_DECIDE,size_active_schur,size_active_schur);CHKERRQ(ierr);
776       ierr = MatSetType(S_Ej_tilda_all,MATAIJ);CHKERRQ(ierr);
777       ierr = MatSeqAIJSetPreallocation(S_Ej_tilda_all,0,nnz);CHKERRQ(ierr);
778       ierr = MatCreate(PETSC_COMM_SELF,&S_Ej_inv_all);CHKERRQ(ierr);
779       ierr = MatSetSizes(S_Ej_inv_all,PETSC_DECIDE,PETSC_DECIDE,size_active_schur,size_active_schur);CHKERRQ(ierr);
780       ierr = MatSetType(S_Ej_inv_all,MATAIJ);CHKERRQ(ierr);
781       ierr = MatSeqAIJSetPreallocation(S_Ej_inv_all,0,nnz);CHKERRQ(ierr);
782     }
783 
784     /* S_Ej_all */
785     cum = cum2 = 0;
786     ierr = MatDenseGetArray(S_all,&S_data);CHKERRQ(ierr);
787     for (i=0;i<sub_schurs->n_subs;i++) {
788       PetscInt j;
789 
790       /* get S_E */
791       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
792       if (mumps_S && sub_schurs->is_hermitian) { /* When using MUMPS data I need to expand to upper triangular (column oriented) */
793         PetscInt k;
794         for (k=0;k<subset_size;k++) {
795           for (j=k;j<subset_size;j++) {
796             work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
797             work[j*subset_size+k] = S_data[cum2+k*size_schur+j];
798           }
799         }
800       } else { /* copy to workspace */
801         PetscInt k;
802         for (k=0;k<subset_size;k++) {
803           for (j=0;j<subset_size;j++) {
804             work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
805           }
806         }
807       }
808       /* insert S_E values */
809       for (j=0;j<subset_size;j++) dummy_idx[j] = cum+j;
810       ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
811 
812       /* if adaptivity is requested, invert S_E block */
813       if (compute_Stilda && ((PetscBTLookup(sub_schurs->is_edge,i) && use_edges) || (!PetscBTLookup(sub_schurs->is_edge,i) && use_faces))) {
814         ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr);
815         ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
816         if (!sub_schurs->is_hermitian) { /* TODO add sytrf/i for symmetric non hermitian */
817           PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,work,&B_N,pivots,&B_ierr));
818           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
819           PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,work,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
820           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
821         } else {
822           PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,work,&B_N,&B_ierr));
823           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
824           PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,work,&B_N,&B_ierr));
825           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
826         }
827         ierr = PetscFPTrapPop();CHKERRQ(ierr);
828         ierr = MatSetValues(S_Ej_inv_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
829       }
830       cum += subset_size;
831       cum2 += subset_size*(size_schur + 1);
832     }
833     ierr = MatDenseRestoreArray(S_all,&S_data);CHKERRQ(ierr);
834 
835 #if defined(PETSC_HAVE_MUMPS)
836     if (mumps_S) {
837       ierr = MatMumpsRestoreSchurComplement(F,&S_all);CHKERRQ(ierr);
838     }
839 #endif
840 
841     if (compute_Stilda) {
842       if (sub_schurs->n_subs == 1 && size_schur == size_active_schur) { /* we already computed the inverse */
843         PetscInt j;
844         for (j=0;j<size_schur;j++) dummy_idx[j] = j;
845         ierr = MatSetValues(S_Ej_tilda_all,size_schur,dummy_idx,size_schur,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
846         ierr = PetscBTSet(sub_schurs->computed_Stilda_subs,0);CHKERRQ(ierr);
847       } else {
848         if (mumps_S) { /* use MatMumps calls to invert S */
849 #if defined(PETSC_HAVE_MUMPS)
850           ierr = MatMumpsInvertSchurComplement(F);CHKERRQ(ierr);
851           ierr = MatMumpsGetSchurComplement(F,&S_all);CHKERRQ(ierr);
852 #endif
853         } else { /* we need to invert explicitly since we are not using MUMPS for S */
854           ierr = MatDenseGetArray(S_all,&S_data);CHKERRQ(ierr);
855           ierr = PetscBLASIntCast(size_schur,&B_N);CHKERRQ(ierr);
856           ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
857           if (!sub_schurs->is_hermitian) { /* TODO add sytrf/i for symmetric non hermitian */
858             PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,S_data,&B_N,pivots,&B_ierr));
859             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
860             PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,S_data,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
861             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
862           } else {
863             PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,S_data,&B_N,&B_ierr));
864             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
865             PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,S_data,&B_N,&B_ierr));
866             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
867           }
868           ierr = PetscFPTrapPop();CHKERRQ(ierr);
869           ierr = MatDenseRestoreArray(S_all,&S_data);CHKERRQ(ierr);
870         }
871         /* S_Ej_tilda_all */
872         cum = cum2 = 0;
873         ierr = MatDenseGetArray(S_all,&S_data);CHKERRQ(ierr);
874         for (i=0;i<sub_schurs->n_subs;i++) {
875           ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
876           /* if adaptivity is requested, invert S_E and insert St_E^-1 blocks */
877           if (((PetscBTLookup(sub_schurs->is_edge,i) && use_edges) || (!PetscBTLookup(sub_schurs->is_edge,i) && use_faces))) {
878             PetscInt j;
879             /* get (St^-1)_E */
880             if (sub_schurs->is_hermitian) { /* I need to expand to upper triangular (column oriented) */
881               PetscInt k;
882               for (k=0;k<subset_size;k++) {
883                 for (j=k;j<subset_size;j++) {
884                   work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
885                   work[j*subset_size+k] = S_data[cum2+k*size_schur+j];
886                 }
887               }
888             } else { /* copy to workspace */
889               PetscInt k;
890               for (k=0;k<subset_size;k++) {
891                 for (j=0;j<subset_size;j++) {
892                   work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
893                 }
894               }
895             }
896             for (j=0;j<subset_size;j++) dummy_idx[j] = cum+j;
897             ierr = MatSetValues(S_Ej_tilda_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
898             ierr = PetscBTSet(sub_schurs->computed_Stilda_subs,i);CHKERRQ(ierr);
899           }
900           cum += subset_size;
901           cum2 += subset_size*(size_schur + 1);
902         }
903         ierr = MatDenseRestoreArray(S_all,&S_data);CHKERRQ(ierr);
904 #if defined(PETSC_HAVE_MUMPS)
905         if (mumps_S) {
906           ierr = MatMumpsRestoreSchurComplement(F,&S_all);CHKERRQ(ierr);
907         }
908 #endif
909       }
910     }
911     ierr = PetscFree2(dummy_idx,work);CHKERRQ(ierr);
912   }
913   ierr = MatDestroy(&F);CHKERRQ(ierr);
914   ierr = ISDestroy(&is_I_layer);CHKERRQ(ierr);
915   ierr = PetscFree(nnz);CHKERRQ(ierr);
916   ierr = MatDestroy(&S_all);CHKERRQ(ierr);
917   ierr = MatDestroy(&A_BB);CHKERRQ(ierr);
918   ierr = MatDestroy(&A_IB);CHKERRQ(ierr);
919   ierr = MatDestroy(&A_BI);CHKERRQ(ierr);
920   ierr = MatAssemblyBegin(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
921   ierr = MatAssemblyEnd(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
922   if (compute_Stilda) {
923     ierr = MatAssemblyBegin(S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
924     ierr = MatAssemblyEnd(S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
925     ierr = MatAssemblyBegin(S_Ej_inv_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
926     ierr = MatAssemblyEnd(S_Ej_inv_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
927   }
928 
929   /* Global matrix of all assembled Schur on subsets */
930   ierr = MatISSetLocalMat(work_mat,sub_schurs->S_Ej_all);CHKERRQ(ierr);
931   ierr = MatISSetMPIXAIJPreallocation_Private(work_mat,global_schur_subsets,PETSC_TRUE);CHKERRQ(ierr);
932   ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
933 
934   /* Get local part of (\sum_j S_Ej) */
935   ierr = ISCreateGeneral(comm_n,local_size,all_local_idx_G,PETSC_OWN_POINTER,&temp_is);CHKERRQ(ierr);
936   ierr = MatDestroy(&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
937   ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
938 
939   /* Compute explicitly (\sum_j S_Ej)^-1 (faster scaling during PCApply, needs extra work when doing setup) */
940   if (faster_deluxe) {
941     Mat         tmpmat;
942     PetscScalar *array;
943     PetscInt    cum;
944 
945     ierr = MatSeqAIJGetArray(sub_schurs->sum_S_Ej_all,&array);CHKERRQ(ierr);
946     cum = 0;
947     for (i=0;i<sub_schurs->n_subs;i++) {
948       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
949       ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr);
950       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
951       if (!sub_schurs->is_hermitian) {
952         PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,array+cum,&B_N,pivots,&B_ierr));
953         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
954         PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,array+cum,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
955         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
956       } else {
957         PetscInt j,k;
958 
959         PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,array+cum,&B_N,&B_ierr));
960         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
961         PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,array+cum,&B_N,&B_ierr));
962         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
963         for (j=0;j<B_N;j++) {
964           for (k=j+1;k<B_N;k++) {
965             array[k*B_N+j+cum] = array[j*B_N+k+cum];
966           }
967         }
968       }
969       ierr = PetscFPTrapPop();CHKERRQ(ierr);
970       cum += subset_size*subset_size;
971     }
972     ierr = MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_all,&array);CHKERRQ(ierr);
973     ierr = MatMatMult(sub_schurs->S_Ej_all,sub_schurs->sum_S_Ej_all,MAT_INITIAL_MATRIX,1.0,&tmpmat);CHKERRQ(ierr);
974     ierr = MatDestroy(&sub_schurs->S_Ej_all);CHKERRQ(ierr);
975     ierr = MatDestroy(&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
976     sub_schurs->S_Ej_all = tmpmat;
977   }
978 
979   /* Get local part of (\sum_j S^-1_Ej) (\sum_j St^-1_Ej) */
980   if (compute_Stilda) {
981     ierr = MatISSetLocalMat(work_mat,S_Ej_tilda_all);CHKERRQ(ierr);
982     ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
983     ierr = MatDestroy(&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
984     ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
985     ierr = MatISSetLocalMat(work_mat,S_Ej_inv_all);CHKERRQ(ierr);
986     ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
987     ierr = MatDestroy(&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
988     ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
989   }
990 
991   /* free workspace */
992   ierr = PetscFree2(Bwork,pivots);CHKERRQ(ierr);
993   ierr = MatDestroy(&global_schur_subsets);CHKERRQ(ierr);
994   ierr = MatDestroy(&S_Ej_tilda_all);CHKERRQ(ierr);
995   ierr = MatDestroy(&S_Ej_inv_all);CHKERRQ(ierr);
996   ierr = MatDestroy(&work_mat);CHKERRQ(ierr);
997   ierr = ISDestroy(&temp_is);CHKERRQ(ierr);
998   ierr = PetscCommDestroy(&comm_n);CHKERRQ(ierr);
999   PetscFunctionReturn(0);
1000 }
1001 
1002 #undef __FUNCT__
1003 #define __FUNCT__ "PCBDDCSubSchursInit"
1004 PetscErrorCode PCBDDCSubSchursInit(PCBDDCSubSchurs sub_schurs, IS is_I, IS is_B, PCBDDCGraph graph, ISLocalToGlobalMapping BtoNmap)
1005 {
1006   IS              *faces,*edges,*all_cc,vertices;
1007   PetscInt        i,n_faces,n_edges,n_all_cc;
1008   PetscBool       is_sorted;
1009   PetscErrorCode  ierr;
1010 
1011   PetscFunctionBegin;
1012   ierr = ISSorted(is_I,&is_sorted);CHKERRQ(ierr);
1013   if (!is_sorted) {
1014     SETERRQ(PetscObjectComm((PetscObject)is_I),PETSC_ERR_PLIB,"IS for I dofs should be shorted");
1015   }
1016   ierr = ISSorted(is_B,&is_sorted);CHKERRQ(ierr);
1017   if (!is_sorted) {
1018     SETERRQ(PetscObjectComm((PetscObject)is_B),PETSC_ERR_PLIB,"IS for B dofs should be shorted");
1019   }
1020 
1021   /* reset any previous data */
1022   ierr = PCBDDCSubSchursReset(sub_schurs);CHKERRQ(ierr);
1023 
1024   /* get index sets for faces and edges (already sorted by global ordering) */
1025   ierr = PCBDDCGraphGetCandidatesIS(graph,&n_faces,&faces,&n_edges,&edges,&vertices);CHKERRQ(ierr);
1026   n_all_cc = n_faces+n_edges;
1027   ierr = PetscBTCreate(n_all_cc,&sub_schurs->is_edge);CHKERRQ(ierr);
1028   ierr = PetscBTCreate(n_all_cc,&sub_schurs->computed_Stilda_subs);CHKERRQ(ierr);
1029   ierr = PetscMalloc1(n_all_cc,&all_cc);CHKERRQ(ierr);
1030   for (i=0;i<n_faces;i++) {
1031     all_cc[i] = faces[i];
1032   }
1033   for (i=0;i<n_edges;i++) {
1034     all_cc[n_faces+i] = edges[i];
1035     ierr = PetscBTSet(sub_schurs->is_edge,n_faces+i);CHKERRQ(ierr);
1036   }
1037   ierr = PetscFree(faces);CHKERRQ(ierr);
1038   ierr = PetscFree(edges);CHKERRQ(ierr);
1039   sub_schurs->is_dir = NULL;
1040   ierr = PCBDDCGraphGetDirichletDofsB(graph,&sub_schurs->is_dir);CHKERRQ(ierr);
1041 
1042   /* Determine if MUMPS can be used */
1043   sub_schurs->use_mumps = PETSC_FALSE;
1044 #if defined(PETSC_HAVE_MUMPS)
1045   sub_schurs->use_mumps = PETSC_TRUE;
1046 #endif
1047 
1048   ierr = PetscObjectReference((PetscObject)is_I);CHKERRQ(ierr);
1049   sub_schurs->is_I = is_I;
1050   ierr = PetscObjectReference((PetscObject)is_B);CHKERRQ(ierr);
1051   sub_schurs->is_B = is_B;
1052   ierr = PetscObjectReference((PetscObject)graph->l2gmap);CHKERRQ(ierr);
1053   sub_schurs->l2gmap = graph->l2gmap;
1054   ierr = PetscObjectReference((PetscObject)BtoNmap);CHKERRQ(ierr);
1055   sub_schurs->BtoNmap = BtoNmap;
1056   sub_schurs->n_subs = n_all_cc;
1057   sub_schurs->is_subs = all_cc;
1058   if (!sub_schurs->use_mumps) { /* sort by local ordering mumps is not present */
1059     for (i=0;i<sub_schurs->n_subs;i++) {
1060       ierr = ISSort(sub_schurs->is_subs[i]);CHKERRQ(ierr);
1061     }
1062   }
1063   sub_schurs->is_vertices = vertices;
1064   sub_schurs->S_Ej_all = NULL;
1065   sub_schurs->sum_S_Ej_all = NULL;
1066   sub_schurs->sum_S_Ej_inv_all = NULL;
1067   sub_schurs->sum_S_Ej_tilda_all = NULL;
1068   sub_schurs->is_Ej_all = NULL;
1069   PetscFunctionReturn(0);
1070 }
1071 
1072 #undef __FUNCT__
1073 #define __FUNCT__ "PCBDDCSubSchursCreate"
1074 PetscErrorCode PCBDDCSubSchursCreate(PCBDDCSubSchurs *sub_schurs)
1075 {
1076   PCBDDCSubSchurs schurs_ctx;
1077   PetscErrorCode  ierr;
1078 
1079   PetscFunctionBegin;
1080   ierr = PetscNew(&schurs_ctx);CHKERRQ(ierr);
1081   schurs_ctx->n_subs = 0;
1082   *sub_schurs = schurs_ctx;
1083   PetscFunctionReturn(0);
1084 }
1085 
1086 #undef __FUNCT__
1087 #define __FUNCT__ "PCBDDCSubSchursDestroy"
1088 PetscErrorCode PCBDDCSubSchursDestroy(PCBDDCSubSchurs *sub_schurs)
1089 {
1090   PetscErrorCode ierr;
1091 
1092   PetscFunctionBegin;
1093   ierr = PCBDDCSubSchursReset(*sub_schurs);CHKERRQ(ierr);
1094   ierr = PetscFree(*sub_schurs);CHKERRQ(ierr);
1095   PetscFunctionReturn(0);
1096 }
1097 
1098 #undef __FUNCT__
1099 #define __FUNCT__ "PCBDDCSubSchursReset"
1100 PetscErrorCode PCBDDCSubSchursReset(PCBDDCSubSchurs sub_schurs)
1101 {
1102   PetscInt       i;
1103   PetscErrorCode ierr;
1104 
1105   PetscFunctionBegin;
1106   ierr = MatDestroy(&sub_schurs->A);CHKERRQ(ierr);
1107   ierr = MatDestroy(&sub_schurs->S);CHKERRQ(ierr);
1108   ierr = ISDestroy(&sub_schurs->is_I);CHKERRQ(ierr);
1109   ierr = ISDestroy(&sub_schurs->is_B);CHKERRQ(ierr);
1110   ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->l2gmap);CHKERRQ(ierr);
1111   ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->BtoNmap);CHKERRQ(ierr);
1112   ierr = MatDestroy(&sub_schurs->S_Ej_all);CHKERRQ(ierr);
1113   ierr = MatDestroy(&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
1114   ierr = MatDestroy(&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
1115   ierr = MatDestroy(&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
1116   ierr = ISDestroy(&sub_schurs->is_Ej_all);CHKERRQ(ierr);
1117   ierr = ISDestroy(&sub_schurs->is_vertices);CHKERRQ(ierr);
1118   ierr = ISDestroy(&sub_schurs->is_dir);CHKERRQ(ierr);
1119   ierr = PetscBTDestroy(&sub_schurs->is_edge);CHKERRQ(ierr);
1120   ierr = PetscBTDestroy(&sub_schurs->computed_Stilda_subs);CHKERRQ(ierr);
1121   for (i=0;i<sub_schurs->n_subs;i++) {
1122     ierr = ISDestroy(&sub_schurs->is_subs[i]);CHKERRQ(ierr);
1123   }
1124   if (sub_schurs->n_subs) {
1125     ierr = PetscFree(sub_schurs->is_subs);CHKERRQ(ierr);
1126   }
1127   if (sub_schurs->reuse_mumps) {
1128     ierr = PCBDDCReuseMumpsReset(sub_schurs->reuse_mumps);CHKERRQ(ierr);
1129   }
1130   ierr = PetscFree(sub_schurs->reuse_mumps);CHKERRQ(ierr);
1131   sub_schurs->n_subs = 0;
1132   PetscFunctionReturn(0);
1133 }
1134 
1135 #undef __FUNCT__
1136 #define __FUNCT__ "PCBDDCAdjGetNextLayer_Private"
1137 PETSC_STATIC_INLINE PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt* queue_tip,PetscInt n_prev,PetscBT touched,PetscInt* xadj,PetscInt* adjncy,PetscInt* n_added)
1138 {
1139   PetscInt       i,j,n;
1140   PetscErrorCode ierr;
1141 
1142   PetscFunctionBegin;
1143   n = 0;
1144   for (i=-n_prev;i<0;i++) {
1145     PetscInt start_dof = queue_tip[i];
1146     for (j=xadj[start_dof];j<xadj[start_dof+1];j++) {
1147       PetscInt dof = adjncy[j];
1148       if (!PetscBTLookup(touched,dof)) {
1149         ierr = PetscBTSet(touched,dof);CHKERRQ(ierr);
1150         queue_tip[n] = dof;
1151         n++;
1152       }
1153     }
1154   }
1155   *n_added = n;
1156   PetscFunctionReturn(0);
1157 }
1158