xref: /petsc/src/ksp/pc/impls/bddc/bddcschurs.c (revision 249c8a7109f3f02db1bdad39df014d76aeb045cb)
1 #include <../src/ksp/pc/impls/bddc/bddc.h>
2 #include <../src/ksp/pc/impls/bddc/bddcprivate.h>
3 #include <../src/mat/impls/dense/seq/dense.h>
4 #include <petscblaslapack.h>
5 
6 PETSC_STATIC_INLINE PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt*,PetscInt,PetscBT,PetscInt*,PetscInt*,PetscInt*);
7 static PetscErrorCode PCBDDCComputeExplicitSchur(Mat,PetscBool,MatReuse,Mat*);
8 static PetscErrorCode PCBDDCReuseSolvers_Interior(PC,Vec,Vec);
9 static PetscErrorCode PCBDDCReuseSolvers_Correction(PC,Vec,Vec);
10 
11 /* if v2 is not present, correction is done in-place */
12 PetscErrorCode PCBDDCReuseSolversBenignAdapt(PCBDDCReuseSolvers ctx, Vec v, Vec v2, PetscBool sol, PetscBool full)
13 {
14   PetscScalar    *array;
15   PetscScalar    *array2;
16   PetscErrorCode ierr;
17 
18   PetscFunctionBegin;
19   if (!ctx->benign_n) PetscFunctionReturn(0);
20   if (sol && full) {
21     PetscInt n_I,size_schur;
22 
23     /* get sizes */
24     ierr = MatGetSize(ctx->benign_csAIB,&size_schur,NULL);CHKERRQ(ierr);
25     ierr = VecGetSize(v,&n_I);CHKERRQ(ierr);
26     n_I = n_I - size_schur;
27     /* get schur sol from array */
28     ierr = VecGetArray(v,&array);CHKERRQ(ierr);
29     ierr = VecPlaceArray(ctx->benign_dummy_schur_vec,array+n_I);CHKERRQ(ierr);
30     ierr = VecRestoreArray(v,&array);CHKERRQ(ierr);
31     /* apply interior sol correction */
32     ierr = MatMultTranspose(ctx->benign_csAIB,ctx->benign_dummy_schur_vec,ctx->benign_corr_work);CHKERRQ(ierr);
33     ierr = VecResetArray(ctx->benign_dummy_schur_vec);CHKERRQ(ierr);
34     ierr = MatMultAdd(ctx->benign_AIIm1ones,ctx->benign_corr_work,v,v);CHKERRQ(ierr);
35   }
36   if (v2) {
37     PetscInt nl;
38 
39     ierr = VecGetArrayRead(v,(const PetscScalar**)&array);CHKERRQ(ierr);
40     ierr = VecGetLocalSize(v2,&nl);CHKERRQ(ierr);
41     ierr = VecGetArray(v2,&array2);CHKERRQ(ierr);
42     ierr = PetscMemcpy(array2,array,nl*sizeof(PetscScalar));CHKERRQ(ierr);
43   } else {
44     ierr = VecGetArray(v,&array);CHKERRQ(ierr);
45     array2 = array;
46   }
47   if (!sol) { /* change rhs */
48     PetscInt n;
49     for (n=0;n<ctx->benign_n;n++) {
50       PetscScalar    sum = 0.;
51       const PetscInt *cols;
52       PetscInt       nz,i;
53 
54       ierr = ISGetLocalSize(ctx->benign_zerodiag_subs[n],&nz);CHKERRQ(ierr);
55       ierr = ISGetIndices(ctx->benign_zerodiag_subs[n],&cols);CHKERRQ(ierr);
56       for (i=0;i<nz-1;i++) sum += array[cols[i]];
57 #if defined(PETSC_USE_COMPLEX)
58       sum = -(PetscRealPart(sum)/nz + PETSC_i*(PetscImaginaryPart(sum)/nz));
59 #else
60       sum = -sum/nz;
61 #endif
62       for (i=0;i<nz-1;i++) array2[cols[i]] += sum;
63       ctx->benign_save_vals[n] = array2[cols[nz-1]];
64       array2[cols[nz-1]] = sum;
65       ierr = ISRestoreIndices(ctx->benign_zerodiag_subs[n],&cols);CHKERRQ(ierr);
66     }
67   } else {
68     PetscInt n;
69     for (n=0;n<ctx->benign_n;n++) {
70       PetscScalar    sum = 0.;
71       const PetscInt *cols;
72       PetscInt       nz,i;
73       ierr = ISGetLocalSize(ctx->benign_zerodiag_subs[n],&nz);CHKERRQ(ierr);
74       ierr = ISGetIndices(ctx->benign_zerodiag_subs[n],&cols);CHKERRQ(ierr);
75       for (i=0;i<nz-1;i++) sum += array[cols[i]];
76 #if defined(PETSC_USE_COMPLEX)
77       sum = -(PetscRealPart(sum)/nz + PETSC_i*(PetscImaginaryPart(sum)/nz));
78 #else
79       sum = -sum/nz;
80 #endif
81       for (i=0;i<nz-1;i++) array2[cols[i]] += sum;
82       array2[cols[nz-1]] = ctx->benign_save_vals[n];
83       ierr = ISRestoreIndices(ctx->benign_zerodiag_subs[n],&cols);CHKERRQ(ierr);
84     }
85   }
86   if (v2) {
87     ierr = VecRestoreArrayRead(v,(const PetscScalar**)&array);CHKERRQ(ierr);
88     ierr = VecRestoreArray(v2,&array2);CHKERRQ(ierr);
89   } else {
90     ierr = VecRestoreArray(v,&array);CHKERRQ(ierr);
91   }
92   if (!sol && full) {
93     Vec      usedv;
94     PetscInt n_I,size_schur;
95 
96     /* get sizes */
97     ierr = MatGetSize(ctx->benign_csAIB,&size_schur,NULL);CHKERRQ(ierr);
98     ierr = VecGetSize(v,&n_I);CHKERRQ(ierr);
99     n_I = n_I - size_schur;
100     /* compute schur rhs correction */
101     if (v2) {
102       usedv = v2;
103     } else {
104       usedv = v;
105     }
106     /* apply schur rhs correction */
107     ierr = MatMultTranspose(ctx->benign_AIIm1ones,usedv,ctx->benign_corr_work);CHKERRQ(ierr);
108     ierr = VecGetArrayRead(usedv,(const PetscScalar**)&array);CHKERRQ(ierr);
109     ierr = VecPlaceArray(ctx->benign_dummy_schur_vec,array+n_I);CHKERRQ(ierr);
110     ierr = VecRestoreArrayRead(usedv,(const PetscScalar**)&array);CHKERRQ(ierr);
111     ierr = MatMultAdd(ctx->benign_csAIB,ctx->benign_corr_work,ctx->benign_dummy_schur_vec,ctx->benign_dummy_schur_vec);CHKERRQ(ierr);
112     ierr = VecResetArray(ctx->benign_dummy_schur_vec);CHKERRQ(ierr);
113   }
114   PetscFunctionReturn(0);
115 }
116 
117 static PetscErrorCode PCBDDCReuseSolvers_Solve_Private(PC pc, Vec rhs, Vec sol, PetscBool transpose, PetscBool full)
118 {
119   PCBDDCReuseSolvers ctx;
120   PetscBool          copy = PETSC_FALSE;
121   PetscErrorCode     ierr;
122 
123   PetscFunctionBegin;
124   ierr = PCShellGetContext(pc,(void **)&ctx);CHKERRQ(ierr);
125   if (full) {
126 #if defined(PETSC_HAVE_MUMPS)
127     ierr = MatMumpsSetIcntl(ctx->F,26,-1);CHKERRQ(ierr);
128 #endif
129 #if defined(PETSC_HAVE_MKL_PARDISO)
130     ierr = MatMkl_PardisoSetCntl(ctx->F,70,0);CHKERRQ(ierr);
131 #endif
132     copy = ctx->has_vertices;
133   } else { /* interior solver */
134 #if defined(PETSC_HAVE_MUMPS)
135     ierr = MatMumpsSetIcntl(ctx->F,26,0);CHKERRQ(ierr);
136 #endif
137 #if defined(PETSC_HAVE_MKL_PARDISO)
138     ierr = MatMkl_PardisoSetCntl(ctx->F,70,1);CHKERRQ(ierr);
139 #endif
140     copy = PETSC_TRUE;
141   }
142   /* copy rhs into factored matrix workspace */
143   if (copy) {
144     PetscInt    n;
145     PetscScalar *array,*array_solver;
146 
147     ierr = VecGetLocalSize(rhs,&n);CHKERRQ(ierr);
148     ierr = VecGetArrayRead(rhs,(const PetscScalar**)&array);CHKERRQ(ierr);
149     ierr = VecGetArray(ctx->rhs,&array_solver);CHKERRQ(ierr);
150     ierr = PetscMemcpy(array_solver,array,n*sizeof(PetscScalar));CHKERRQ(ierr);
151     ierr = VecRestoreArray(ctx->rhs,&array_solver);CHKERRQ(ierr);
152     ierr = VecRestoreArrayRead(rhs,(const PetscScalar**)&array);CHKERRQ(ierr);
153 
154     ierr = PCBDDCReuseSolversBenignAdapt(ctx,ctx->rhs,NULL,PETSC_FALSE,full);CHKERRQ(ierr);
155     if (transpose) {
156       ierr = MatSolveTranspose(ctx->F,ctx->rhs,ctx->sol);CHKERRQ(ierr);
157     } else {
158       ierr = MatSolve(ctx->F,ctx->rhs,ctx->sol);CHKERRQ(ierr);
159     }
160     ierr = PCBDDCReuseSolversBenignAdapt(ctx,ctx->sol,NULL,PETSC_TRUE,full);CHKERRQ(ierr);
161 
162     /* get back data to caller worskpace */
163     ierr = VecGetArrayRead(ctx->sol,(const PetscScalar**)&array_solver);CHKERRQ(ierr);
164     ierr = VecGetArray(sol,&array);CHKERRQ(ierr);
165     ierr = PetscMemcpy(array,array_solver,n*sizeof(PetscScalar));CHKERRQ(ierr);
166     ierr = VecRestoreArray(sol,&array);CHKERRQ(ierr);
167     ierr = VecRestoreArrayRead(ctx->sol,(const PetscScalar**)&array_solver);CHKERRQ(ierr);
168   } else {
169     if (ctx->benign_n) {
170       ierr = PCBDDCReuseSolversBenignAdapt(ctx,rhs,ctx->rhs,PETSC_FALSE,full);CHKERRQ(ierr);
171       if (transpose) {
172         ierr = MatSolveTranspose(ctx->F,ctx->rhs,sol);CHKERRQ(ierr);
173       } else {
174         ierr = MatSolve(ctx->F,ctx->rhs,sol);CHKERRQ(ierr);
175       }
176       ierr = PCBDDCReuseSolversBenignAdapt(ctx,sol,NULL,PETSC_TRUE,full);CHKERRQ(ierr);
177     } else {
178       if (transpose) {
179         ierr = MatSolveTranspose(ctx->F,rhs,sol);CHKERRQ(ierr);
180       } else {
181         ierr = MatSolve(ctx->F,rhs,sol);CHKERRQ(ierr);
182       }
183     }
184   }
185   /* restore defaults */
186 #if defined(PETSC_HAVE_MUMPS)
187   ierr = MatMumpsSetIcntl(ctx->F,26,-1);CHKERRQ(ierr);
188 #endif
189 #if defined(PETSC_HAVE_MKL_PARDISO)
190   ierr = MatMkl_PardisoSetCntl(ctx->F,70,0);CHKERRQ(ierr);
191 #endif
192   PetscFunctionReturn(0);
193 }
194 
195 static PetscErrorCode PCBDDCReuseSolvers_Correction(PC pc, Vec rhs, Vec sol)
196 {
197   PetscErrorCode   ierr;
198 
199   PetscFunctionBegin;
200   ierr = PCBDDCReuseSolvers_Solve_Private(pc,rhs,sol,PETSC_FALSE,PETSC_TRUE);CHKERRQ(ierr);
201   PetscFunctionReturn(0);
202 }
203 
204 static PetscErrorCode PCBDDCReuseSolvers_CorrectionTranspose(PC pc, Vec rhs, Vec sol)
205 {
206   PetscErrorCode   ierr;
207 
208   PetscFunctionBegin;
209   ierr = PCBDDCReuseSolvers_Solve_Private(pc,rhs,sol,PETSC_TRUE,PETSC_TRUE);CHKERRQ(ierr);
210   PetscFunctionReturn(0);
211 }
212 
213 static PetscErrorCode PCBDDCReuseSolvers_Interior(PC pc, Vec rhs, Vec sol)
214 {
215   PetscErrorCode   ierr;
216 
217   PetscFunctionBegin;
218   ierr = PCBDDCReuseSolvers_Solve_Private(pc,rhs,sol,PETSC_FALSE,PETSC_FALSE);CHKERRQ(ierr);
219   PetscFunctionReturn(0);
220 }
221 
222 static PetscErrorCode PCBDDCReuseSolvers_InteriorTranspose(PC pc, Vec rhs, Vec sol)
223 {
224   PetscErrorCode   ierr;
225 
226   PetscFunctionBegin;
227   ierr = PCBDDCReuseSolvers_Solve_Private(pc,rhs,sol,PETSC_TRUE,PETSC_FALSE);CHKERRQ(ierr);
228   PetscFunctionReturn(0);
229 }
230 
231 static PetscErrorCode PCBDDCReuseSolvers_View(PC pc, PetscViewer viewer)
232 {
233   PCBDDCReuseSolvers ctx;
234   PetscBool          iascii;
235   PetscErrorCode     ierr;
236 
237   PetscFunctionBegin;
238   ierr = PCShellGetContext(pc,(void **)&ctx);CHKERRQ(ierr);
239   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
240   if (iascii) {
241     ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr);
242   }
243   ierr = MatView(ctx->F,viewer);CHKERRQ(ierr);
244   if (iascii) {
245     ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
246   }
247   PetscFunctionReturn(0);
248 }
249 
250 static PetscErrorCode PCBDDCReuseSolversReset(PCBDDCReuseSolvers reuse)
251 {
252   PetscInt       i;
253   PetscErrorCode ierr;
254 
255   PetscFunctionBegin;
256   ierr = MatDestroy(&reuse->F);CHKERRQ(ierr);
257   ierr = VecDestroy(&reuse->sol);CHKERRQ(ierr);
258   ierr = VecDestroy(&reuse->rhs);CHKERRQ(ierr);
259   ierr = PCDestroy(&reuse->interior_solver);CHKERRQ(ierr);
260   ierr = PCDestroy(&reuse->correction_solver);CHKERRQ(ierr);
261   ierr = ISDestroy(&reuse->is_R);CHKERRQ(ierr);
262   ierr = ISDestroy(&reuse->is_B);CHKERRQ(ierr);
263   ierr = VecScatterDestroy(&reuse->correction_scatter_B);CHKERRQ(ierr);
264   ierr = VecDestroy(&reuse->sol_B);CHKERRQ(ierr);
265   ierr = VecDestroy(&reuse->rhs_B);CHKERRQ(ierr);
266   for (i=0;i<reuse->benign_n;i++) {
267     ierr = ISDestroy(&reuse->benign_zerodiag_subs[i]);CHKERRQ(ierr);
268   }
269   ierr = PetscFree(reuse->benign_zerodiag_subs);CHKERRQ(ierr);
270   ierr = PetscFree(reuse->benign_save_vals);CHKERRQ(ierr);
271   ierr = MatDestroy(&reuse->benign_csAIB);CHKERRQ(ierr);
272   ierr = MatDestroy(&reuse->benign_AIIm1ones);CHKERRQ(ierr);
273   ierr = VecDestroy(&reuse->benign_corr_work);CHKERRQ(ierr);
274   ierr = VecDestroy(&reuse->benign_dummy_schur_vec);CHKERRQ(ierr);
275   PetscFunctionReturn(0);
276 }
277 
278 static PetscErrorCode PCBDDCComputeExplicitSchur(Mat M, PetscBool issym, MatReuse reuse, Mat *S)
279 {
280   Mat            B, C, D, Bd, Cd, AinvBd;
281   KSP            ksp;
282   PC             pc;
283   PetscBool      isLU, isILU, isCHOL, Bdense, Cdense;
284   PetscReal      fill = 2.0;
285   PetscInt       n_I;
286   PetscMPIInt    size;
287   PetscErrorCode ierr;
288 
289   PetscFunctionBegin;
290   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)M),&size);CHKERRQ(ierr);
291   if (size != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not for parallel matrices");
292   if (reuse == MAT_REUSE_MATRIX) {
293     PetscBool Sdense;
294 
295     ierr = PetscObjectTypeCompare((PetscObject)*S, MATSEQDENSE, &Sdense);CHKERRQ(ierr);
296     if (!Sdense) SETERRQ(PetscObjectComm((PetscObject)M),PETSC_ERR_SUP,"S should dense");
297   }
298   ierr = MatSchurComplementGetSubMatrices(M, NULL, NULL, &B, &C, &D);CHKERRQ(ierr);
299   ierr = MatSchurComplementGetKSP(M, &ksp);CHKERRQ(ierr);
300   ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr);
301   ierr = PetscObjectTypeCompare((PetscObject) pc, PCLU, &isLU);CHKERRQ(ierr);
302   ierr = PetscObjectTypeCompare((PetscObject) pc, PCILU, &isILU);CHKERRQ(ierr);
303   ierr = PetscObjectTypeCompare((PetscObject) pc, PCCHOLESKY, &isCHOL);CHKERRQ(ierr);
304   ierr = PetscObjectTypeCompare((PetscObject) B, MATSEQDENSE, &Bdense);CHKERRQ(ierr);
305   ierr = PetscObjectTypeCompare((PetscObject) C, MATSEQDENSE, &Cdense);CHKERRQ(ierr);
306   ierr = MatGetSize(B,&n_I,NULL);CHKERRQ(ierr);
307   if (n_I) {
308     if (!Bdense) {
309       ierr = MatConvert(B, MATSEQDENSE, MAT_INITIAL_MATRIX, &Bd);CHKERRQ(ierr);
310     } else {
311       Bd = B;
312     }
313 
314     if (isLU || isILU || isCHOL) {
315       Mat fact;
316       ierr = KSPSetUp(ksp);CHKERRQ(ierr);
317       ierr = PCFactorGetMatrix(pc, &fact);CHKERRQ(ierr);
318       ierr = MatDuplicate(Bd, MAT_DO_NOT_COPY_VALUES, &AinvBd);CHKERRQ(ierr);
319       ierr = MatMatSolve(fact, Bd, AinvBd);CHKERRQ(ierr);
320     } else {
321       PetscBool ex = PETSC_TRUE;
322 
323       if (ex) {
324         Mat Ainvd;
325 
326         ierr = PCComputeExplicitOperator(pc, &Ainvd);CHKERRQ(ierr);
327         ierr = MatMatMult(Ainvd, Bd, MAT_INITIAL_MATRIX, fill, &AinvBd);CHKERRQ(ierr);
328         ierr = MatDestroy(&Ainvd);CHKERRQ(ierr);
329       } else {
330         Vec         sol,rhs;
331         PetscScalar *arrayrhs,*arraysol;
332         PetscInt    i,nrhs,n;
333 
334         ierr = MatDuplicate(Bd, MAT_DO_NOT_COPY_VALUES, &AinvBd);CHKERRQ(ierr);
335         ierr = MatGetSize(Bd,&n,&nrhs);CHKERRQ(ierr);
336         ierr = MatDenseGetArray(Bd,&arrayrhs);CHKERRQ(ierr);
337         ierr = MatDenseGetArray(AinvBd,&arraysol);CHKERRQ(ierr);
338         ierr = KSPGetSolution(ksp,&sol);CHKERRQ(ierr);
339         ierr = KSPGetRhs(ksp,&rhs);CHKERRQ(ierr);
340         for (i=0;i<nrhs;i++) {
341           ierr = VecPlaceArray(rhs,arrayrhs+i*n);CHKERRQ(ierr);
342           ierr = VecPlaceArray(sol,arraysol+i*n);CHKERRQ(ierr);
343           ierr = KSPSolve(ksp,rhs,sol);CHKERRQ(ierr);
344           ierr = VecResetArray(rhs);CHKERRQ(ierr);
345           ierr = VecResetArray(sol);CHKERRQ(ierr);
346         }
347         ierr = MatDenseRestoreArray(Bd,&arrayrhs);CHKERRQ(ierr);
348         ierr = MatDenseRestoreArray(AinvBd,&arrayrhs);CHKERRQ(ierr);
349       }
350     }
351     if (!Bdense & !issym) {
352       ierr = MatDestroy(&Bd);CHKERRQ(ierr);
353     }
354 
355     if (!issym) {
356       if (!Cdense) {
357         ierr = MatConvert(C, MATSEQDENSE, MAT_INITIAL_MATRIX, &Cd);CHKERRQ(ierr);
358       } else {
359         Cd = C;
360       }
361       ierr = MatMatMult(Cd, AinvBd, reuse, fill, S);CHKERRQ(ierr);
362       if (!Cdense) {
363         ierr = MatDestroy(&Cd);CHKERRQ(ierr);
364       }
365     } else {
366       ierr = MatTransposeMatMult(Bd, AinvBd, reuse, fill, S);CHKERRQ(ierr);
367       if (!Bdense) {
368         ierr = MatDestroy(&Bd);CHKERRQ(ierr);
369       }
370     }
371     ierr = MatDestroy(&AinvBd);CHKERRQ(ierr);
372   }
373 
374   if (D) {
375     Mat       Dd;
376     PetscBool Ddense;
377 
378     ierr = PetscObjectTypeCompare((PetscObject)D,MATSEQDENSE,&Ddense);CHKERRQ(ierr);
379     if (!Ddense) {
380       ierr = MatConvert(D, MATSEQDENSE, MAT_INITIAL_MATRIX, &Dd);CHKERRQ(ierr);
381     } else {
382       Dd = D;
383     }
384     if (n_I) {
385       ierr = MatAYPX(*S,-1.0,Dd,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
386     } else {
387       if (reuse == MAT_INITIAL_MATRIX) {
388         ierr = MatDuplicate(Dd,MAT_COPY_VALUES,S);CHKERRQ(ierr);
389       } else {
390         ierr = MatCopy(Dd,*S,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
391       }
392     }
393     if (!Ddense) {
394       ierr = MatDestroy(&Dd);CHKERRQ(ierr);
395     }
396   } else {
397     ierr = MatScale(*S,-1.0);CHKERRQ(ierr);
398   }
399   PetscFunctionReturn(0);
400 }
401 
402 PetscErrorCode PCBDDCSubSchursSetUp(PCBDDCSubSchurs sub_schurs, Mat Ain, Mat Sin, PetscBool exact_schur, PetscInt xadj[], PetscInt adjncy[], PetscInt nlayers, Vec scaling, PetscBool compute_Stilda, PetscBool reuse_solvers, PetscBool benign_trick, PetscInt benign_n, PetscInt benign_p0_lidx[], IS benign_zerodiag_subs[], Mat change, IS change_primal)
403 {
404   Mat                    F,A_II,A_IB,A_BI,A_BB,AE_II;
405   Mat                    S_all;
406   Mat                    global_schur_subsets,work_mat,*submats;
407   ISLocalToGlobalMapping l2gmap_subsets;
408   IS                     is_I,is_I_layer;
409   IS                     all_subsets,all_subsets_mult,all_subsets_n;
410   PetscInt               *nnz,*all_local_idx_N;
411   PetscInt               *auxnum1,*auxnum2;
412   PetscInt               i,subset_size,max_subset_size;
413   PetscInt               n_B,extra,local_size,global_size;
414   PetscBLASInt           B_N,B_ierr,B_lwork,*pivots;
415   PetscScalar            *Bwork;
416   MPI_Comm               comm_n;
417   PetscBool              deluxe = PETSC_TRUE;
418   PetscBool              use_potr = PETSC_FALSE, use_sytr = PETSC_FALSE;
419   PetscBool              print_schurs = PETSC_FALSE;
420   PetscErrorCode         ierr;
421 
422   PetscFunctionBegin;
423   ierr = MatDestroy(&sub_schurs->A);CHKERRQ(ierr);
424   ierr = MatDestroy(&sub_schurs->S);CHKERRQ(ierr);
425   /* convert matrix if needed */
426   if (Ain) {
427     PetscBool isseqaij;
428     ierr = PetscObjectTypeCompare((PetscObject)Ain,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
429     if (isseqaij) {
430       ierr = PetscObjectReference((PetscObject)Ain);CHKERRQ(ierr);
431       sub_schurs->A = Ain;
432     } else {
433       ierr = MatConvert(Ain,MATSEQAIJ,MAT_INITIAL_MATRIX,&sub_schurs->A);CHKERRQ(ierr);
434     }
435   }
436 
437   ierr = PetscObjectReference((PetscObject)Sin);CHKERRQ(ierr);
438   sub_schurs->S = Sin;
439   if (sub_schurs->schur_explicit) {
440     sub_schurs->schur_explicit = (PetscBool)(!!sub_schurs->A);
441   }
442 
443   /* preliminary checks */
444   if (!sub_schurs->schur_explicit && compute_Stilda) SETERRQ(PetscObjectComm((PetscObject)sub_schurs->l2gmap),PETSC_ERR_SUP,"Adaptive selection of constraints requires MUMPS and/or MKL_PARDISO");
445 
446   if (benign_trick) sub_schurs->is_posdef = PETSC_FALSE;
447 
448   /* debug prints */
449   if (sub_schurs->debug) {
450     PetscMPIInt size,rank;
451     PetscInt    nr,*print_schurs_ranks;
452     PetscBool   flg;
453 
454     ierr = MPI_Comm_size(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&size);CHKERRQ(ierr);
455     ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&rank);CHKERRQ(ierr);
456     nr   = size;
457     ierr = PetscMalloc1(nr,&print_schurs_ranks);CHKERRQ(ierr);
458     ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)sub_schurs->l2gmap),sub_schurs->prefix,"BDDC sub_schurs options","PC");CHKERRQ(ierr);
459     ierr = PetscOptionsIntArray("-sub_schurs_debug_ranks","Ranks to debug (all if the option is not used)",NULL,print_schurs_ranks,&nr,&flg);CHKERRQ(ierr);
460     if (!flg) print_schurs = PETSC_TRUE;
461     else {
462       for (i=0;i<nr;i++) if (print_schurs_ranks[i] == (PetscInt)rank) { print_schurs = PETSC_TRUE; break; }
463     }
464     ierr = PetscOptionsEnd();CHKERRQ(ierr);
465     ierr = PetscFree(print_schurs_ranks);CHKERRQ(ierr);
466   }
467 
468 
469   /* restrict work on active processes */
470   if (sub_schurs->restrict_comm) {
471     PetscSubcomm subcomm;
472     PetscMPIInt  color,rank;
473 
474     color = 0;
475     if (!sub_schurs->n_subs) color = 1; /* this can happen if we are in a multilevel case or if the subdomain is disconnected */
476     ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&rank);CHKERRQ(ierr);
477     ierr = PetscSubcommCreate(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&subcomm);CHKERRQ(ierr);
478     ierr = PetscSubcommSetNumber(subcomm,2);CHKERRQ(ierr);
479     ierr = PetscSubcommSetTypeGeneral(subcomm,color,rank);CHKERRQ(ierr);
480     ierr = PetscCommDuplicate(PetscSubcommChild(subcomm),&comm_n,NULL);CHKERRQ(ierr);
481     ierr = PetscSubcommDestroy(&subcomm);CHKERRQ(ierr);
482     if (!sub_schurs->n_subs) {
483       ierr = PetscCommDestroy(&comm_n);CHKERRQ(ierr);
484       PetscFunctionReturn(0);
485     }
486   } else {
487     ierr = PetscCommDuplicate(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&comm_n,NULL);CHKERRQ(ierr);
488   }
489 
490   /* get Schur complement matrices */
491   if (!sub_schurs->schur_explicit) {
492     Mat       tA_IB,tA_BI,tA_BB;
493     PetscBool isseqsbaij;
494     ierr = MatSchurComplementGetSubMatrices(sub_schurs->S,&A_II,NULL,&tA_IB,&tA_BI,&tA_BB);CHKERRQ(ierr);
495     ierr = PetscObjectTypeCompare((PetscObject)tA_BB,MATSEQSBAIJ,&isseqsbaij);CHKERRQ(ierr);
496     if (isseqsbaij) {
497       ierr = MatConvert(tA_BB,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_BB);CHKERRQ(ierr);
498       ierr = MatConvert(tA_IB,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_IB);CHKERRQ(ierr);
499       ierr = MatConvert(tA_BI,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_BI);CHKERRQ(ierr);
500     } else {
501       ierr = PetscObjectReference((PetscObject)tA_BB);CHKERRQ(ierr);
502       A_BB = tA_BB;
503       ierr = PetscObjectReference((PetscObject)tA_IB);CHKERRQ(ierr);
504       A_IB = tA_IB;
505       ierr = PetscObjectReference((PetscObject)tA_BI);CHKERRQ(ierr);
506       A_BI = tA_BI;
507     }
508   } else {
509     A_II = NULL;
510     A_IB = NULL;
511     A_BI = NULL;
512     A_BB = NULL;
513   }
514   S_all = NULL;
515 
516   /* determine interior problems */
517   ierr = ISGetLocalSize(sub_schurs->is_I,&i);CHKERRQ(ierr);
518   if (nlayers >= 0 && i) { /* Interior problems can be different from the original one */
519     PetscBT                touched;
520     const PetscInt*        idx_B;
521     PetscInt               n_I,n_B,n_local_dofs,n_prev_added,j,layer,*local_numbering;
522 
523     if (!xadj) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cannot request layering without adjacency");
524     /* get sizes */
525     ierr = ISGetLocalSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr);
526     ierr = ISGetLocalSize(sub_schurs->is_B,&n_B);CHKERRQ(ierr);
527 
528     ierr = PetscMalloc1(n_I+n_B,&local_numbering);CHKERRQ(ierr);
529     ierr = PetscBTCreate(n_I+n_B,&touched);CHKERRQ(ierr);
530     ierr = PetscBTMemzero(n_I+n_B,touched);CHKERRQ(ierr);
531 
532     /* all boundary dofs must be skipped when adding layers */
533     ierr = ISGetIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr);
534     for (j=0;j<n_B;j++) {
535       ierr = PetscBTSet(touched,idx_B[j]);CHKERRQ(ierr);
536     }
537     ierr = PetscMemcpy(local_numbering,idx_B,n_B*sizeof(PetscInt));CHKERRQ(ierr);
538     ierr = ISRestoreIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr);
539 
540     /* add prescribed number of layers of dofs */
541     n_local_dofs = n_B;
542     n_prev_added = n_B;
543     for (layer=0;layer<nlayers;layer++) {
544       PetscInt n_added;
545       if (n_local_dofs == n_I+n_B) break;
546       if (n_local_dofs > n_I+n_B) 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);
547       ierr = PCBDDCAdjGetNextLayer_Private(local_numbering+n_local_dofs,n_prev_added,touched,xadj,adjncy,&n_added);CHKERRQ(ierr);
548       n_prev_added = n_added;
549       n_local_dofs += n_added;
550       if (!n_added) break;
551     }
552     ierr = PetscBTDestroy(&touched);CHKERRQ(ierr);
553 
554     /* IS for I layer dofs in original numbering */
555     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);
556     ierr = PetscFree(local_numbering);CHKERRQ(ierr);
557     ierr = ISSort(is_I_layer);CHKERRQ(ierr);
558     /* IS for I layer dofs in I numbering */
559     if (!sub_schurs->schur_explicit) {
560       ISLocalToGlobalMapping ItoNmap;
561       ierr = ISLocalToGlobalMappingCreateIS(sub_schurs->is_I,&ItoNmap);CHKERRQ(ierr);
562       ierr = ISGlobalToLocalMappingApplyIS(ItoNmap,IS_GTOLM_DROP,is_I_layer,&is_I);CHKERRQ(ierr);
563       ierr = ISLocalToGlobalMappingDestroy(&ItoNmap);CHKERRQ(ierr);
564 
565       /* II block */
566       ierr = MatCreateSubMatrix(A_II,is_I,is_I,MAT_INITIAL_MATRIX,&AE_II);CHKERRQ(ierr);
567     }
568   } else {
569     PetscInt n_I;
570 
571     /* IS for I dofs in original numbering */
572     ierr = PetscObjectReference((PetscObject)sub_schurs->is_I);CHKERRQ(ierr);
573     is_I_layer = sub_schurs->is_I;
574 
575     /* IS for I dofs in I numbering (strided 1) */
576     if (!sub_schurs->schur_explicit) {
577       ierr = ISGetSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr);
578       ierr = ISCreateStride(PetscObjectComm((PetscObject)sub_schurs->is_I),n_I,0,1,&is_I);CHKERRQ(ierr);
579 
580       /* II block is the same */
581       ierr = PetscObjectReference((PetscObject)A_II);CHKERRQ(ierr);
582       AE_II = A_II;
583     }
584   }
585 
586   /* Get info on subset sizes and sum of all subsets sizes */
587   max_subset_size = 0;
588   local_size = 0;
589   for (i=0;i<sub_schurs->n_subs;i++) {
590     ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
591     max_subset_size = PetscMax(subset_size,max_subset_size);
592     local_size += subset_size;
593   }
594 
595   /* Work arrays for local indices */
596   extra = 0;
597   ierr = ISGetLocalSize(sub_schurs->is_B,&n_B);CHKERRQ(ierr);
598   if (sub_schurs->schur_explicit && is_I_layer) {
599     ierr = ISGetLocalSize(is_I_layer,&extra);CHKERRQ(ierr);
600   }
601   ierr = PetscMalloc1(n_B+extra,&all_local_idx_N);CHKERRQ(ierr);
602   if (extra) {
603     const PetscInt *idxs;
604     ierr = ISGetIndices(is_I_layer,&idxs);CHKERRQ(ierr);
605     ierr = PetscMemcpy(all_local_idx_N,idxs,extra*sizeof(PetscInt));CHKERRQ(ierr);
606     ierr = ISRestoreIndices(is_I_layer,&idxs);CHKERRQ(ierr);
607   }
608   ierr = PetscMalloc1(local_size,&nnz);CHKERRQ(ierr);
609   ierr = PetscMalloc1(sub_schurs->n_subs,&auxnum1);CHKERRQ(ierr);
610   ierr = PetscMalloc1(sub_schurs->n_subs,&auxnum2);CHKERRQ(ierr);
611 
612   /* Get local indices in local numbering */
613   local_size = 0;
614   for (i=0;i<sub_schurs->n_subs;i++) {
615     PetscInt j;
616     const    PetscInt *idxs;
617 
618     ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
619     ierr = ISGetIndices(sub_schurs->is_subs[i],&idxs);CHKERRQ(ierr);
620     /* start (smallest in global ordering) and multiplicity */
621     auxnum1[i] = idxs[0];
622     auxnum2[i] = subset_size;
623     /* subset indices in local numbering */
624     ierr = PetscMemcpy(all_local_idx_N+local_size+extra,idxs,subset_size*sizeof(PetscInt));CHKERRQ(ierr);
625     ierr = ISRestoreIndices(sub_schurs->is_subs[i],&idxs);CHKERRQ(ierr);
626     for (j=0;j<subset_size;j++) nnz[local_size+j] = subset_size;
627     local_size += subset_size;
628   }
629 
630   /* allocate extra workspace needed only for GETRI or SYTRF */
631   use_potr = use_sytr = PETSC_FALSE;
632   if (benign_trick || (sub_schurs->is_hermitian && sub_schurs->is_posdef)) {
633     use_potr = PETSC_TRUE;
634   } else if (sub_schurs->is_symmetric) {
635     use_sytr = PETSC_TRUE;
636   }
637   if (local_size && !use_potr) {
638     PetscScalar  lwork,dummyscalar = 0.;
639     PetscBLASInt dummyint = 0;
640 
641     B_lwork = -1;
642     ierr = PetscBLASIntCast(local_size,&B_N);CHKERRQ(ierr);
643     ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
644     if (use_sytr) {
645       PetscStackCallBLAS("LAPACKsytrf",LAPACKsytrf_("L",&B_N,&dummyscalar,&B_N,&dummyint,&lwork,&B_lwork,&B_ierr));
646       if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to SYTRF Lapack routine %d",(int)B_ierr);
647     } else {
648       PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,&dummyscalar,&B_N,&dummyint,&lwork,&B_lwork,&B_ierr));
649       if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to GETRI Lapack routine %d",(int)B_ierr);
650     }
651     ierr = PetscFPTrapPop();CHKERRQ(ierr);
652     ierr = PetscBLASIntCast((PetscInt)PetscRealPart(lwork),&B_lwork);CHKERRQ(ierr);
653     ierr = PetscMalloc2(B_lwork,&Bwork,B_N,&pivots);CHKERRQ(ierr);
654   } else {
655     Bwork = NULL;
656     pivots = NULL;
657   }
658 
659   /* prepare parallel matrices for summing up properly schurs on subsets */
660   ierr = ISCreateGeneral(comm_n,sub_schurs->n_subs,auxnum1,PETSC_OWN_POINTER,&all_subsets_n);CHKERRQ(ierr);
661   ierr = ISLocalToGlobalMappingApplyIS(sub_schurs->l2gmap,all_subsets_n,&all_subsets);CHKERRQ(ierr);
662   ierr = ISDestroy(&all_subsets_n);CHKERRQ(ierr);
663   ierr = ISCreateGeneral(comm_n,sub_schurs->n_subs,auxnum2,PETSC_OWN_POINTER,&all_subsets_mult);CHKERRQ(ierr);
664   ierr = ISRenumber(all_subsets,all_subsets_mult,&global_size,&all_subsets_n);CHKERRQ(ierr);
665   ierr = ISDestroy(&all_subsets);CHKERRQ(ierr);
666   ierr = ISDestroy(&all_subsets_mult);CHKERRQ(ierr);
667   ierr = ISGetLocalSize(all_subsets_n,&i);CHKERRQ(ierr);
668   if (i != local_size) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Invalid size of new subset! %D != %D",i,local_size);
669   ierr = ISLocalToGlobalMappingCreateIS(all_subsets_n,&l2gmap_subsets);CHKERRQ(ierr);
670   ierr = MatCreateIS(comm_n,1,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size,l2gmap_subsets,NULL,&work_mat);CHKERRQ(ierr);
671   ierr = ISLocalToGlobalMappingDestroy(&l2gmap_subsets);CHKERRQ(ierr);
672   ierr = MatCreate(PetscObjectComm((PetscObject)work_mat),&global_schur_subsets);CHKERRQ(ierr);
673   ierr = MatSetSizes(global_schur_subsets,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size);CHKERRQ(ierr);
674   ierr = MatSetType(global_schur_subsets,MATMPIAIJ);CHKERRQ(ierr);
675 
676   /* subset indices in local boundary numbering */
677   if (!sub_schurs->is_Ej_all) {
678     PetscInt *all_local_idx_B;
679 
680     ierr = PetscMalloc1(local_size,&all_local_idx_B);CHKERRQ(ierr);
681     ierr = ISGlobalToLocalMappingApply(sub_schurs->BtoNmap,IS_GTOLM_DROP,local_size,all_local_idx_N+extra,&subset_size,all_local_idx_B);CHKERRQ(ierr);
682     if (subset_size != local_size) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in sub_schurs serial (BtoNmap)! %D != %D",subset_size,local_size);
683     ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size,all_local_idx_B,PETSC_OWN_POINTER,&sub_schurs->is_Ej_all);CHKERRQ(ierr);
684   }
685 
686   if (change) {
687     ISLocalToGlobalMapping BtoS;
688     IS                     change_primal_B;
689     IS                     change_primal_all;
690 
691     if (sub_schurs->change_primal_sub) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"This should not happen");
692     if (sub_schurs->change) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"This should not happen");
693     ierr = PetscMalloc1(sub_schurs->n_subs,&sub_schurs->change_primal_sub);CHKERRQ(ierr);
694     for (i=0;i<sub_schurs->n_subs;i++) {
695       ISLocalToGlobalMapping NtoS;
696       ierr = ISLocalToGlobalMappingCreateIS(sub_schurs->is_subs[i],&NtoS);CHKERRQ(ierr);
697       ierr = ISGlobalToLocalMappingApplyIS(NtoS,IS_GTOLM_DROP,change_primal,&sub_schurs->change_primal_sub[i]);CHKERRQ(ierr);
698       ierr = ISLocalToGlobalMappingDestroy(&NtoS);CHKERRQ(ierr);
699     }
700     ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,change_primal,&change_primal_B);CHKERRQ(ierr);
701     ierr = ISLocalToGlobalMappingCreateIS(sub_schurs->is_Ej_all,&BtoS);CHKERRQ(ierr);
702     ierr = ISGlobalToLocalMappingApplyIS(BtoS,IS_GTOLM_DROP,change_primal_B,&change_primal_all);CHKERRQ(ierr);
703     ierr = ISLocalToGlobalMappingDestroy(&BtoS);CHKERRQ(ierr);
704     ierr = ISDestroy(&change_primal_B);CHKERRQ(ierr);
705     ierr = PetscMalloc1(sub_schurs->n_subs,&sub_schurs->change);CHKERRQ(ierr);
706     for (i=0;i<sub_schurs->n_subs;i++) {
707       Mat change_sub;
708 
709       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
710       ierr = KSPCreate(PETSC_COMM_SELF,&sub_schurs->change[i]);CHKERRQ(ierr);
711       ierr = KSPSetType(sub_schurs->change[i],KSPPREONLY);CHKERRQ(ierr);
712       if (!sub_schurs->change_with_qr) {
713         ierr = MatCreateSubMatrix(change,sub_schurs->is_subs[i],sub_schurs->is_subs[i],MAT_INITIAL_MATRIX,&change_sub);CHKERRQ(ierr);
714       } else {
715         Mat change_subt;
716         ierr = MatCreateSubMatrix(change,sub_schurs->is_subs[i],sub_schurs->is_subs[i],MAT_INITIAL_MATRIX,&change_subt);CHKERRQ(ierr);
717         ierr = MatConvert(change_subt,MATSEQDENSE,MAT_INITIAL_MATRIX,&change_sub);CHKERRQ(ierr);
718         ierr = MatDestroy(&change_subt);CHKERRQ(ierr);
719       }
720       ierr = KSPSetOperators(sub_schurs->change[i],change_sub,change_sub);CHKERRQ(ierr);
721       ierr = MatDestroy(&change_sub);CHKERRQ(ierr);
722       ierr = KSPSetOptionsPrefix(sub_schurs->change[i],sub_schurs->prefix);CHKERRQ(ierr);
723       ierr = KSPAppendOptionsPrefix(sub_schurs->change[i],"sub_schurs_change_");CHKERRQ(ierr);
724     }
725     ierr = ISDestroy(&change_primal_all);CHKERRQ(ierr);
726   }
727 
728   /* Local matrix of all local Schur on subsets (transposed) */
729   if (!sub_schurs->S_Ej_all) {
730     ierr = MatCreate(PETSC_COMM_SELF,&sub_schurs->S_Ej_all);CHKERRQ(ierr);
731     ierr = MatSetSizes(sub_schurs->S_Ej_all,PETSC_DECIDE,PETSC_DECIDE,local_size,local_size);CHKERRQ(ierr);
732     ierr = MatSetType(sub_schurs->S_Ej_all,MATAIJ);CHKERRQ(ierr);
733     ierr = MatSeqAIJSetPreallocation(sub_schurs->S_Ej_all,0,nnz);CHKERRQ(ierr);
734   }
735 
736   /* Compute Schur complements explicitly */
737   F = NULL;
738   if (!sub_schurs->schur_explicit) {
739     /* this code branch is used when MatFactor with Schur complement support is not present or when explicitly requested;
740        it is not efficient, unless the economic version of the scaling is used */
741     Mat         S_Ej_expl;
742     PetscScalar *work;
743     PetscInt    j,*dummy_idx;
744     PetscBool   Sdense;
745 
746     ierr = PetscMalloc2(max_subset_size,&dummy_idx,max_subset_size*max_subset_size,&work);CHKERRQ(ierr);
747     local_size = 0;
748     for (i=0;i<sub_schurs->n_subs;i++) {
749       IS  is_subset_B;
750       Mat AE_EE,AE_IE,AE_EI,S_Ej;
751 
752       /* subsets in original and boundary numbering */
753       ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,sub_schurs->is_subs[i],&is_subset_B);CHKERRQ(ierr);
754       /* EE block */
755       ierr = MatCreateSubMatrix(A_BB,is_subset_B,is_subset_B,MAT_INITIAL_MATRIX,&AE_EE);CHKERRQ(ierr);
756       /* IE block */
757       ierr = MatCreateSubMatrix(A_IB,is_I,is_subset_B,MAT_INITIAL_MATRIX,&AE_IE);CHKERRQ(ierr);
758       /* EI block */
759       if (sub_schurs->is_symmetric) {
760         ierr = MatCreateTranspose(AE_IE,&AE_EI);CHKERRQ(ierr);
761       } else if (sub_schurs->is_hermitian) {
762         ierr = MatCreateHermitianTranspose(AE_IE,&AE_EI);CHKERRQ(ierr);
763       } else {
764         ierr = MatCreateSubMatrix(A_BI,is_subset_B,is_I,MAT_INITIAL_MATRIX,&AE_EI);CHKERRQ(ierr);
765       }
766       ierr = ISDestroy(&is_subset_B);CHKERRQ(ierr);
767       ierr = MatCreateSchurComplement(AE_II,AE_II,AE_IE,AE_EI,AE_EE,&S_Ej);CHKERRQ(ierr);
768       ierr = MatDestroy(&AE_EE);CHKERRQ(ierr);
769       ierr = MatDestroy(&AE_IE);CHKERRQ(ierr);
770       ierr = MatDestroy(&AE_EI);CHKERRQ(ierr);
771       if (AE_II == A_II) { /* we can reuse the same ksp */
772         KSP ksp;
773         ierr = MatSchurComplementGetKSP(sub_schurs->S,&ksp);CHKERRQ(ierr);
774         ierr = MatSchurComplementSetKSP(S_Ej,ksp);CHKERRQ(ierr);
775       } else { /* build new ksp object which inherits ksp and pc types from the original one */
776         KSP       origksp,schurksp;
777         PC        origpc,schurpc;
778         KSPType   ksp_type;
779         PetscInt  n_internal;
780         PetscBool ispcnone;
781 
782         ierr = MatSchurComplementGetKSP(sub_schurs->S,&origksp);CHKERRQ(ierr);
783         ierr = MatSchurComplementGetKSP(S_Ej,&schurksp);CHKERRQ(ierr);
784         ierr = KSPGetType(origksp,&ksp_type);CHKERRQ(ierr);
785         ierr = KSPSetType(schurksp,ksp_type);CHKERRQ(ierr);
786         ierr = KSPGetPC(schurksp,&schurpc);CHKERRQ(ierr);
787         ierr = KSPGetPC(origksp,&origpc);CHKERRQ(ierr);
788         ierr = PetscObjectTypeCompare((PetscObject)origpc,PCNONE,&ispcnone);CHKERRQ(ierr);
789         if (!ispcnone) {
790           PCType pc_type;
791           ierr = PCGetType(origpc,&pc_type);CHKERRQ(ierr);
792           ierr = PCSetType(schurpc,pc_type);CHKERRQ(ierr);
793         } else {
794           ierr = PCSetType(schurpc,PCLU);CHKERRQ(ierr);
795         }
796         ierr = ISGetSize(is_I,&n_internal);CHKERRQ(ierr);
797         if (n_internal) { /* UMFPACK gives error with 0 sized problems */
798           MatSolverType solver = NULL;
799           ierr = PCFactorGetMatSolverType(origpc,(MatSolverType*)&solver);CHKERRQ(ierr);
800           if (solver) {
801             ierr = PCFactorSetMatSolverType(schurpc,solver);CHKERRQ(ierr);
802           }
803         }
804         ierr = KSPSetUp(schurksp);CHKERRQ(ierr);
805       }
806       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
807       ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&S_Ej_expl);CHKERRQ(ierr);
808       ierr = PCBDDCComputeExplicitSchur(S_Ej,sub_schurs->is_symmetric,MAT_REUSE_MATRIX,&S_Ej_expl);CHKERRQ(ierr);
809       ierr = PetscObjectTypeCompare((PetscObject)S_Ej_expl,MATSEQDENSE,&Sdense);CHKERRQ(ierr);
810       if (Sdense) {
811         for (j=0;j<subset_size;j++) {
812           dummy_idx[j]=local_size+j;
813         }
814         ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
815       } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not yet implemented for sparse matrices");
816       ierr = MatDestroy(&S_Ej);CHKERRQ(ierr);
817       ierr = MatDestroy(&S_Ej_expl);CHKERRQ(ierr);
818       local_size += subset_size;
819     }
820     ierr = PetscFree2(dummy_idx,work);CHKERRQ(ierr);
821     /* free */
822     ierr = ISDestroy(&is_I);CHKERRQ(ierr);
823     ierr = MatDestroy(&AE_II);CHKERRQ(ierr);
824     ierr = PetscFree(all_local_idx_N);CHKERRQ(ierr);
825   } else {
826     Mat         A,cs_AIB_mat = NULL,benign_AIIm1_ones_mat = NULL;
827     Vec         Dall = NULL;
828     IS          is_A_all,*is_p_r = NULL;
829     PetscScalar *work,*S_data,*schur_factor,infty = PETSC_MAX_REAL;
830     PetscInt    n,n_I,*dummy_idx,size_schur,size_active_schur,cum,cum2;
831     PetscBool   economic,solver_S,S_lower_triangular = PETSC_FALSE;
832     PetscBool   schur_has_vertices,factor_workaround;
833     PetscBool   use_cholesky;
834 
835     /* get sizes */
836     n_I = 0;
837     if (is_I_layer) {
838       ierr = ISGetLocalSize(is_I_layer,&n_I);CHKERRQ(ierr);
839     }
840     economic = PETSC_FALSE;
841     ierr = ISGetLocalSize(sub_schurs->is_I,&cum);CHKERRQ(ierr);
842     if (cum != n_I) economic = PETSC_TRUE;
843     ierr = MatGetLocalSize(sub_schurs->A,&n,NULL);CHKERRQ(ierr);
844     size_active_schur = local_size;
845 
846     /* import scaling vector (wrong formulation if we have 3D edges) */
847     if (scaling && compute_Stilda) {
848       const PetscScalar *array;
849       PetscScalar       *array2;
850       const PetscInt    *idxs;
851       PetscInt          i;
852 
853       ierr = ISGetIndices(sub_schurs->is_Ej_all,&idxs);CHKERRQ(ierr);
854       ierr = VecCreateSeq(PETSC_COMM_SELF,size_active_schur,&Dall);CHKERRQ(ierr);
855       ierr = VecGetArrayRead(scaling,&array);CHKERRQ(ierr);
856       ierr = VecGetArray(Dall,&array2);CHKERRQ(ierr);
857       for (i=0;i<size_active_schur;i++) array2[i] = array[idxs[i]];
858       ierr = VecRestoreArray(Dall,&array2);CHKERRQ(ierr);
859       ierr = VecRestoreArrayRead(scaling,&array);CHKERRQ(ierr);
860       ierr = ISRestoreIndices(sub_schurs->is_Ej_all,&idxs);CHKERRQ(ierr);
861       deluxe = PETSC_FALSE;
862     }
863 
864     /* size active schurs does not count any dirichlet or vertex dof on the interface */
865     factor_workaround = PETSC_FALSE;
866     schur_has_vertices = PETSC_FALSE;
867     cum = n_I+size_active_schur;
868     if (sub_schurs->is_dir) {
869       const PetscInt* idxs;
870       PetscInt        n_dir;
871 
872       ierr = ISGetLocalSize(sub_schurs->is_dir,&n_dir);CHKERRQ(ierr);
873       ierr = ISGetIndices(sub_schurs->is_dir,&idxs);CHKERRQ(ierr);
874       ierr = PetscMemcpy(all_local_idx_N+cum,idxs,n_dir*sizeof(PetscInt));CHKERRQ(ierr);
875       ierr = ISRestoreIndices(sub_schurs->is_dir,&idxs);CHKERRQ(ierr);
876       cum += n_dir;
877       factor_workaround = PETSC_TRUE;
878     }
879     /* include the primal vertices in the Schur complement */
880     if (exact_schur && sub_schurs->is_vertices && (compute_Stilda || benign_n)) {
881       PetscInt n_v;
882 
883       ierr = ISGetLocalSize(sub_schurs->is_vertices,&n_v);CHKERRQ(ierr);
884       if (n_v) {
885         const PetscInt* idxs;
886 
887         ierr = ISGetIndices(sub_schurs->is_vertices,&idxs);CHKERRQ(ierr);
888         ierr = PetscMemcpy(all_local_idx_N+cum,idxs,n_v*sizeof(PetscInt));CHKERRQ(ierr);
889         ierr = ISRestoreIndices(sub_schurs->is_vertices,&idxs);CHKERRQ(ierr);
890         cum += n_v;
891         factor_workaround = PETSC_TRUE;
892         schur_has_vertices = PETSC_TRUE;
893       }
894     }
895     size_schur = cum - n_I;
896     ierr = ISCreateGeneral(PETSC_COMM_SELF,cum,all_local_idx_N,PETSC_OWN_POINTER,&is_A_all);CHKERRQ(ierr);
897     if (cum == n) {
898       ierr = ISSetPermutation(is_A_all);CHKERRQ(ierr);
899       ierr = MatPermute(sub_schurs->A,is_A_all,is_A_all,&A);CHKERRQ(ierr);
900     } else {
901       ierr = MatCreateSubMatrix(sub_schurs->A,is_A_all,is_A_all,MAT_INITIAL_MATRIX,&A);CHKERRQ(ierr);
902     }
903     ierr = MatSetOptionsPrefix(A,sub_schurs->prefix);CHKERRQ(ierr);
904     ierr = MatAppendOptionsPrefix(A,"sub_schurs_");CHKERRQ(ierr);
905 
906     /* if we actually change the basis for the pressures, LDL^T factors will use a lot of memory
907        n^2 instead of n^1.5 or something. This is a workaround */
908     if (benign_n) {
909       Vec                    v;
910       ISLocalToGlobalMapping N_to_reor;
911       IS                     is_p0,is_p0_p;
912       PetscScalar            *cs_AIB,*AIIm1_data;
913       PetscInt               sizeA;
914 
915       ierr = ISLocalToGlobalMappingCreateIS(is_A_all,&N_to_reor);CHKERRQ(ierr);
916       ierr = ISCreateGeneral(PETSC_COMM_SELF,benign_n,benign_p0_lidx,PETSC_COPY_VALUES,&is_p0);CHKERRQ(ierr);
917       ierr = ISGlobalToLocalMappingApplyIS(N_to_reor,IS_GTOLM_DROP,is_p0,&is_p0_p);CHKERRQ(ierr);
918       ierr = ISDestroy(&is_p0);CHKERRQ(ierr);
919       ierr = MatCreateVecs(A,&v,NULL);CHKERRQ(ierr);
920       ierr = VecGetSize(v,&sizeA);CHKERRQ(ierr);
921       ierr = MatCreateSeqDense(PETSC_COMM_SELF,sizeA,benign_n,NULL,&benign_AIIm1_ones_mat);CHKERRQ(ierr);
922       ierr = MatCreateSeqDense(PETSC_COMM_SELF,size_schur,benign_n,NULL,&cs_AIB_mat);CHKERRQ(ierr);
923       ierr = MatDenseGetArray(cs_AIB_mat,&cs_AIB);CHKERRQ(ierr);
924       ierr = MatDenseGetArray(benign_AIIm1_ones_mat,&AIIm1_data);CHKERRQ(ierr);
925       ierr = PetscMalloc1(benign_n,&is_p_r);CHKERRQ(ierr);
926       /* compute colsum of A_IB restricted to pressures */
927       for (i=0;i<benign_n;i++) {
928         Vec            benign_AIIm1_ones;
929         PetscScalar    *array;
930         const PetscInt *idxs;
931         PetscInt       j,nz;
932 
933         ierr = ISGlobalToLocalMappingApplyIS(N_to_reor,IS_GTOLM_DROP,benign_zerodiag_subs[i],&is_p_r[i]);CHKERRQ(ierr);
934         ierr = ISGetLocalSize(is_p_r[i],&nz);CHKERRQ(ierr);
935         ierr = ISGetIndices(is_p_r[i],&idxs);CHKERRQ(ierr);
936         for (j=0;j<nz;j++) AIIm1_data[idxs[j]+sizeA*i] = 1.;
937         ierr = ISRestoreIndices(is_p_r[i],&idxs);CHKERRQ(ierr);
938         ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,sizeA,AIIm1_data+sizeA*i,&benign_AIIm1_ones);CHKERRQ(ierr);
939         ierr = MatMult(A,benign_AIIm1_ones,v);CHKERRQ(ierr);
940         ierr = VecDestroy(&benign_AIIm1_ones);CHKERRQ(ierr);
941         ierr = VecGetArray(v,&array);CHKERRQ(ierr);
942         for (j=0;j<size_schur;j++) {
943 #if defined(PETSC_USE_COMPLEX)
944           cs_AIB[i*size_schur + j] = (PetscRealPart(array[j+n_I])/nz + PETSC_i*(PetscImaginaryPart(array[j+n_I])/nz));
945 #else
946           cs_AIB[i*size_schur + j] = array[j+n_I]/nz;
947 #endif
948         }
949         ierr = VecRestoreArray(v,&array);CHKERRQ(ierr);
950       }
951       ierr = MatDenseRestoreArray(cs_AIB_mat,&cs_AIB);CHKERRQ(ierr);
952       ierr = MatDenseRestoreArray(benign_AIIm1_ones_mat,&AIIm1_data);CHKERRQ(ierr);
953       ierr = VecDestroy(&v);CHKERRQ(ierr);
954       ierr = MatSetOption(A,MAT_KEEP_NONZERO_PATTERN,PETSC_FALSE);CHKERRQ(ierr);
955       ierr = MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr);
956       ierr = MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr);
957       ierr = MatZeroRowsColumnsIS(A,is_p0_p,1.0,NULL,NULL);CHKERRQ(ierr);
958       ierr = ISDestroy(&is_p0_p);CHKERRQ(ierr);
959       ierr = ISLocalToGlobalMappingDestroy(&N_to_reor);CHKERRQ(ierr);
960     }
961     ierr = MatSetOption(A,MAT_SYMMETRIC,sub_schurs->is_symmetric);CHKERRQ(ierr);
962     ierr = MatSetOption(A,MAT_HERMITIAN,sub_schurs->is_hermitian);CHKERRQ(ierr);
963     ierr = MatSetOption(A,MAT_SPD,sub_schurs->is_posdef);CHKERRQ(ierr);
964 
965     /* for complexes, symmetric and hermitian at the same time implies null imaginary part */
966     use_cholesky = (PetscBool)((use_potr || use_sytr) && sub_schurs->is_hermitian && sub_schurs->is_symmetric);
967 
968     /* when using the benign subspace trick, the local Schur complements are SPD */
969     if (benign_trick) sub_schurs->is_posdef = PETSC_TRUE;
970 
971     if (n_I) {
972       IS is_schur;
973 
974       if (use_cholesky) {
975         ierr = MatGetFactor(A,sub_schurs->mat_solver_type,MAT_FACTOR_CHOLESKY,&F);CHKERRQ(ierr);
976       } else {
977         ierr = MatGetFactor(A,sub_schurs->mat_solver_type,MAT_FACTOR_LU,&F);CHKERRQ(ierr);
978       }
979       ierr = MatSetErrorIfFailure(A,PETSC_TRUE);CHKERRQ(ierr);
980 
981       /* subsets ordered last */
982       ierr = ISCreateStride(PETSC_COMM_SELF,size_schur,n_I,1,&is_schur);CHKERRQ(ierr);
983       ierr = MatFactorSetSchurIS(F,is_schur);CHKERRQ(ierr);
984       ierr = ISDestroy(&is_schur);CHKERRQ(ierr);
985 
986       /* factorization step */
987       if (use_cholesky) {
988         ierr = MatCholeskyFactorSymbolic(F,A,NULL,NULL);CHKERRQ(ierr);
989 #if defined(PETSC_HAVE_MUMPS) /* be sure that icntl 19 is not set by command line */
990         ierr = MatMumpsSetIcntl(F,19,2);CHKERRQ(ierr);
991 #endif
992         ierr = MatCholeskyFactorNumeric(F,A,NULL);CHKERRQ(ierr);
993         S_lower_triangular = PETSC_TRUE;
994       } else {
995         ierr = MatLUFactorSymbolic(F,A,NULL,NULL,NULL);CHKERRQ(ierr);
996 #if defined(PETSC_HAVE_MUMPS) /* be sure that icntl 19 is not set by command line */
997         ierr = MatMumpsSetIcntl(F,19,3);CHKERRQ(ierr);
998 #endif
999         ierr = MatLUFactorNumeric(F,A,NULL);CHKERRQ(ierr);
1000       }
1001       ierr = MatViewFromOptions(F,(PetscObject)A,"-mat_factor_view");CHKERRQ(ierr);
1002 
1003       if (print_schurs) {
1004         PetscViewer viewer;
1005         char        filename[256];
1006         Mat         S;
1007         IS          is;
1008 
1009         ierr = PetscSNPrintf(filename,sizeof(filename),"sub_schurs_Schur_r%d.m",PetscGlobalRank);CHKERRQ(ierr);
1010         ierr = PetscViewerASCIIOpen(PETSC_COMM_SELF,filename,&viewer);CHKERRQ(ierr);
1011         ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);CHKERRQ(ierr);
1012         ierr = PetscObjectSetName((PetscObject)A,"A");CHKERRQ(ierr);
1013         ierr = MatView(A,viewer);CHKERRQ(ierr);
1014         ierr = MatFactorCreateSchurComplement(F,&S,NULL);CHKERRQ(ierr);
1015         ierr = PetscObjectSetName((PetscObject)S,"S");CHKERRQ(ierr);
1016         ierr = MatView(S,viewer);CHKERRQ(ierr);
1017         ierr = MatDestroy(&S);CHKERRQ(ierr);
1018         ierr = ISCreateStride(PETSC_COMM_SELF,n_I,0,1,&is);CHKERRQ(ierr);
1019         ierr = PetscObjectSetName((PetscObject)is,"I");CHKERRQ(ierr);
1020         ierr = ISView(is,viewer);CHKERRQ(ierr);
1021         ierr = ISDestroy(&is);CHKERRQ(ierr);
1022         ierr = ISCreateStride(PETSC_COMM_SELF,size_schur,n_I,1,&is);CHKERRQ(ierr);
1023         ierr = PetscObjectSetName((PetscObject)is,"B");CHKERRQ(ierr);
1024         ierr = ISView(is,viewer);CHKERRQ(ierr);
1025         ierr = ISDestroy(&is);CHKERRQ(ierr);
1026         ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
1027       }
1028 
1029       /* get explicit Schur Complement computed during numeric factorization */
1030       ierr = MatFactorGetSchurComplement(F,&S_all,NULL);CHKERRQ(ierr);
1031       ierr = MatSetOption(S_all,MAT_SPD,sub_schurs->is_posdef);CHKERRQ(ierr);
1032       ierr = MatSetOption(S_all,MAT_HERMITIAN,sub_schurs->is_hermitian);CHKERRQ(ierr);
1033 
1034       /* we can reuse the solvers if we are not using the economic version */
1035       reuse_solvers = (PetscBool)(reuse_solvers && !economic);
1036       factor_workaround = (PetscBool)(reuse_solvers && factor_workaround);
1037       if (!sub_schurs->is_posdef && factor_workaround && compute_Stilda && size_active_schur)
1038         reuse_solvers = factor_workaround = PETSC_FALSE;
1039 
1040       solver_S = PETSC_TRUE;
1041 
1042       /* update the Schur complement with the change of basis on the pressures */
1043       if (benign_n) {
1044         PetscScalar *S_data,*cs_AIB,*AIIm1_data;
1045         Vec         v;
1046         PetscInt    sizeA;
1047 
1048         ierr = MatDenseGetArray(S_all,&S_data);CHKERRQ(ierr);
1049         ierr = MatCreateVecs(A,&v,NULL);CHKERRQ(ierr);
1050         ierr = VecGetSize(v,&sizeA);CHKERRQ(ierr);
1051 #if defined(PETSC_HAVE_MUMPS)
1052         ierr = MatMumpsSetIcntl(F,26,0);CHKERRQ(ierr);
1053 #endif
1054 #if defined(PETSC_HAVE_MKL_PARDISO)
1055         ierr = MatMkl_PardisoSetCntl(F,70,1);CHKERRQ(ierr);
1056 #endif
1057         ierr = MatDenseGetArray(cs_AIB_mat,&cs_AIB);CHKERRQ(ierr);
1058         ierr = MatDenseGetArray(benign_AIIm1_ones_mat,&AIIm1_data);CHKERRQ(ierr);
1059         for (i=0;i<benign_n;i++) {
1060           Vec            benign_AIIm1_ones;
1061           PetscScalar    *array,sum = 0.,one = 1.;
1062           const PetscInt *idxs;
1063           PetscInt       j,nz;
1064           PetscBLASInt   B_k,B_n;
1065 
1066           ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,sizeA,AIIm1_data+sizeA*i,&benign_AIIm1_ones);CHKERRQ(ierr);
1067           ierr = VecCopy(benign_AIIm1_ones,v);CHKERRQ(ierr);
1068           ierr = MatSolve(F,v,benign_AIIm1_ones);CHKERRQ(ierr);
1069           ierr = ISGetLocalSize(is_p_r[i],&nz);CHKERRQ(ierr);
1070           ierr = ISGetIndices(is_p_r[i],&idxs);CHKERRQ(ierr);
1071           /* p0 dof (eliminated) is excluded from the sum */
1072           for (j=0;j<nz-1;j++) sum -= AIIm1_data[idxs[j]+sizeA*i];
1073           ierr = ISRestoreIndices(is_p_r[i],&idxs);CHKERRQ(ierr);
1074           ierr = MatMult(A,benign_AIIm1_ones,v);CHKERRQ(ierr);
1075           ierr = VecGetArray(v,&array);CHKERRQ(ierr);
1076           /* perform sparse rank updates on symmetric Schur (TODO: move outside of the loop?) */
1077           /* cs_AIB already scaled by 1./nz */
1078           B_k = 1;
1079           ierr = PetscBLASIntCast(size_schur,&B_n);CHKERRQ(ierr);
1080           PetscStackCallBLAS("BLASsyrk",BLASsyrk_("L","N",&B_n,&B_k,&sum,cs_AIB+i*size_schur,&B_n,&one,S_data,&B_n));
1081           sum  = 1.;
1082           PetscStackCallBLAS("BLASsyr2k",BLASsyr2k_("L","N",&B_n,&B_k,&sum,array+n_I,&B_n,cs_AIB+i*size_schur,&B_n,&one,S_data,&B_n));
1083           ierr = VecRestoreArray(v,&array);CHKERRQ(ierr);
1084           /* set p0 entry of AIIm1_ones to zero */
1085           ierr = ISGetIndices(is_p_r[i],&idxs);CHKERRQ(ierr);
1086           for (j=0;j<benign_n;j++) AIIm1_data[idxs[nz-1]+sizeA*j] = 0.;
1087           ierr = ISRestoreIndices(is_p_r[i],&idxs);CHKERRQ(ierr);
1088           ierr = VecDestroy(&benign_AIIm1_ones);CHKERRQ(ierr);
1089         }
1090         if (!S_lower_triangular) { /* I need to expand the upper triangular data (column oriented) */
1091           PetscInt k,j;
1092           for (k=0;k<size_schur;k++) {
1093             for (j=k;j<size_schur;j++) {
1094               S_data[j*size_schur+k] = PetscConj(S_data[k*size_schur+j]);
1095             }
1096           }
1097         }
1098 
1099         /* restore defaults */
1100 #if defined(PETSC_HAVE_MUMPS)
1101         ierr = MatMumpsSetIcntl(F,26,-1);CHKERRQ(ierr);
1102 #endif
1103 #if defined(PETSC_HAVE_MKL_PARDISO)
1104         ierr = MatMkl_PardisoSetCntl(F,70,0);CHKERRQ(ierr);
1105 #endif
1106         ierr = MatDenseRestoreArray(cs_AIB_mat,&cs_AIB);CHKERRQ(ierr);
1107         ierr = MatDenseRestoreArray(benign_AIIm1_ones_mat,&AIIm1_data);CHKERRQ(ierr);
1108         ierr = VecDestroy(&v);CHKERRQ(ierr);
1109         ierr = MatDenseRestoreArray(S_all,&S_data);CHKERRQ(ierr);
1110       }
1111       if (!reuse_solvers) {
1112         for (i=0;i<benign_n;i++) {
1113           ierr = ISDestroy(&is_p_r[i]);CHKERRQ(ierr);
1114         }
1115         ierr = PetscFree(is_p_r);CHKERRQ(ierr);
1116         ierr = MatDestroy(&cs_AIB_mat);CHKERRQ(ierr);
1117         ierr = MatDestroy(&benign_AIIm1_ones_mat);CHKERRQ(ierr);
1118       }
1119     } else { /* we can't use MatFactor when size_schur == size_of_the_problem */
1120       ierr = MatConvert(A,MATSEQDENSE,MAT_INITIAL_MATRIX,&S_all);CHKERRQ(ierr);
1121       reuse_solvers = PETSC_FALSE; /* TODO: why we can't reuse the solvers here? */
1122       factor_workaround = PETSC_FALSE;
1123       solver_S = PETSC_FALSE;
1124     }
1125 
1126     if (reuse_solvers) {
1127       Mat                A_II,Afake;
1128       Vec                vec1_B;
1129       PCBDDCReuseSolvers msolv_ctx;
1130       PetscInt           n_R;
1131 
1132       if (sub_schurs->reuse_solver) {
1133         ierr = PCBDDCReuseSolversReset(sub_schurs->reuse_solver);CHKERRQ(ierr);
1134       } else {
1135         ierr = PetscNew(&sub_schurs->reuse_solver);CHKERRQ(ierr);
1136       }
1137       msolv_ctx = sub_schurs->reuse_solver;
1138       ierr = MatSchurComplementGetSubMatrices(sub_schurs->S,&A_II,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
1139       ierr = PetscObjectReference((PetscObject)F);CHKERRQ(ierr);
1140       msolv_ctx->F = F;
1141       ierr = MatCreateVecs(F,&msolv_ctx->sol,NULL);CHKERRQ(ierr);
1142       /* currently PETSc has no support for MatSolve(F,x,x), so cheat and let rhs and sol share the same memory */
1143       {
1144         PetscScalar *array;
1145         PetscInt    n;
1146 
1147         ierr = VecGetLocalSize(msolv_ctx->sol,&n);CHKERRQ(ierr);
1148         ierr = VecGetArray(msolv_ctx->sol,&array);CHKERRQ(ierr);
1149         ierr = VecCreateSeqWithArray(PetscObjectComm((PetscObject)msolv_ctx->sol),1,n,array,&msolv_ctx->rhs);CHKERRQ(ierr);
1150         ierr = VecRestoreArray(msolv_ctx->sol,&array);CHKERRQ(ierr);
1151       }
1152       msolv_ctx->has_vertices = schur_has_vertices;
1153 
1154       /* interior solver */
1155       ierr = PCCreate(PetscObjectComm((PetscObject)A_II),&msolv_ctx->interior_solver);CHKERRQ(ierr);
1156       ierr = PCSetOperators(msolv_ctx->interior_solver,A_II,A_II);CHKERRQ(ierr);
1157       ierr = PCSetType(msolv_ctx->interior_solver,PCSHELL);CHKERRQ(ierr);
1158       ierr = PCShellSetName(msolv_ctx->interior_solver,"Interior solver (w/o Schur factorization)");CHKERRQ(ierr);
1159       ierr = PCShellSetContext(msolv_ctx->interior_solver,msolv_ctx);CHKERRQ(ierr);
1160       ierr = PCShellSetView(msolv_ctx->interior_solver,PCBDDCReuseSolvers_View);CHKERRQ(ierr);
1161       ierr = PCShellSetApply(msolv_ctx->interior_solver,PCBDDCReuseSolvers_Interior);CHKERRQ(ierr);
1162       ierr = PCShellSetApplyTranspose(msolv_ctx->interior_solver,PCBDDCReuseSolvers_InteriorTranspose);CHKERRQ(ierr);
1163 
1164       /* correction solver */
1165       ierr = PCCreate(PetscObjectComm((PetscObject)A_II),&msolv_ctx->correction_solver);CHKERRQ(ierr);
1166       ierr = PCSetType(msolv_ctx->correction_solver,PCSHELL);CHKERRQ(ierr);
1167       ierr = PCShellSetName(msolv_ctx->correction_solver,"Correction solver (with Schur factorization)");CHKERRQ(ierr);
1168       ierr = PCShellSetContext(msolv_ctx->correction_solver,msolv_ctx);CHKERRQ(ierr);
1169       ierr = PCShellSetView(msolv_ctx->interior_solver,PCBDDCReuseSolvers_View);CHKERRQ(ierr);
1170       ierr = PCShellSetApply(msolv_ctx->correction_solver,PCBDDCReuseSolvers_Correction);CHKERRQ(ierr);
1171       ierr = PCShellSetApplyTranspose(msolv_ctx->correction_solver,PCBDDCReuseSolvers_CorrectionTranspose);CHKERRQ(ierr);
1172 
1173       /* scatter and vecs for Schur complement solver */
1174       ierr = MatCreateVecs(S_all,&msolv_ctx->sol_B,&msolv_ctx->rhs_B);CHKERRQ(ierr);
1175       ierr = MatCreateVecs(sub_schurs->S,&vec1_B,NULL);CHKERRQ(ierr);
1176       if (!schur_has_vertices) {
1177         ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,is_A_all,&msolv_ctx->is_B);CHKERRQ(ierr);
1178         ierr = VecScatterCreate(vec1_B,msolv_ctx->is_B,msolv_ctx->sol_B,NULL,&msolv_ctx->correction_scatter_B);CHKERRQ(ierr);
1179         ierr = PetscObjectReference((PetscObject)is_A_all);CHKERRQ(ierr);
1180         msolv_ctx->is_R = is_A_all;
1181       } else {
1182         IS              is_B_all;
1183         const PetscInt* idxs;
1184         PetscInt        dual,n_v,n;
1185 
1186         ierr = ISGetLocalSize(sub_schurs->is_vertices,&n_v);CHKERRQ(ierr);
1187         dual = size_schur - n_v;
1188         ierr = ISGetLocalSize(is_A_all,&n);CHKERRQ(ierr);
1189         ierr = ISGetIndices(is_A_all,&idxs);CHKERRQ(ierr);
1190         ierr = ISCreateGeneral(PetscObjectComm((PetscObject)is_A_all),dual,idxs+n_I,PETSC_COPY_VALUES,&is_B_all);CHKERRQ(ierr);
1191         ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,is_B_all,&msolv_ctx->is_B);CHKERRQ(ierr);
1192         ierr = ISDestroy(&is_B_all);CHKERRQ(ierr);
1193         ierr = ISCreateStride(PetscObjectComm((PetscObject)is_A_all),dual,0,1,&is_B_all);CHKERRQ(ierr);
1194         ierr = VecScatterCreate(vec1_B,msolv_ctx->is_B,msolv_ctx->sol_B,is_B_all,&msolv_ctx->correction_scatter_B);CHKERRQ(ierr);
1195         ierr = ISDestroy(&is_B_all);CHKERRQ(ierr);
1196         ierr = ISCreateGeneral(PetscObjectComm((PetscObject)is_A_all),n-n_v,idxs,PETSC_COPY_VALUES,&msolv_ctx->is_R);CHKERRQ(ierr);
1197         ierr = ISRestoreIndices(is_A_all,&idxs);CHKERRQ(ierr);
1198       }
1199       ierr = ISGetLocalSize(msolv_ctx->is_R,&n_R);CHKERRQ(ierr);
1200       ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,n_R,n_R,0,NULL,&Afake);CHKERRQ(ierr);
1201       ierr = MatAssemblyBegin(Afake,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1202       ierr = MatAssemblyEnd(Afake,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1203       ierr = PCSetOperators(msolv_ctx->correction_solver,Afake,Afake);CHKERRQ(ierr);
1204       ierr = MatDestroy(&Afake);CHKERRQ(ierr);
1205       ierr = VecDestroy(&vec1_B);CHKERRQ(ierr);
1206 
1207       /* communicate benign info to solver context */
1208       if (benign_n) {
1209         PetscScalar *array;
1210 
1211         msolv_ctx->benign_n = benign_n;
1212         msolv_ctx->benign_zerodiag_subs = is_p_r;
1213         ierr = PetscMalloc1(benign_n,&msolv_ctx->benign_save_vals);CHKERRQ(ierr);
1214         msolv_ctx->benign_csAIB = cs_AIB_mat;
1215         ierr = MatCreateVecs(cs_AIB_mat,&msolv_ctx->benign_corr_work,NULL);CHKERRQ(ierr);
1216         ierr = VecGetArray(msolv_ctx->benign_corr_work,&array);CHKERRQ(ierr);
1217         ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,size_schur,array,&msolv_ctx->benign_dummy_schur_vec);CHKERRQ(ierr);
1218         ierr = VecRestoreArray(msolv_ctx->benign_corr_work,&array);CHKERRQ(ierr);
1219         msolv_ctx->benign_AIIm1ones = benign_AIIm1_ones_mat;
1220       }
1221     } else {
1222       if (sub_schurs->reuse_solver) {
1223         ierr = PCBDDCReuseSolversReset(sub_schurs->reuse_solver);CHKERRQ(ierr);
1224       }
1225       ierr = PetscFree(sub_schurs->reuse_solver);CHKERRQ(ierr);
1226     }
1227     ierr = MatDestroy(&A);CHKERRQ(ierr);
1228     ierr = ISDestroy(&is_A_all);CHKERRQ(ierr);
1229 
1230     /* Work arrays */
1231     ierr = PetscMalloc2(max_subset_size,&dummy_idx,max_subset_size*max_subset_size,&work);CHKERRQ(ierr);
1232 
1233     /* matrices for deluxe scaling and adaptive selection */
1234     if (compute_Stilda) {
1235       if (!sub_schurs->sum_S_Ej_tilda_all) {
1236         ierr = MatDuplicate(sub_schurs->S_Ej_all,MAT_SHARE_NONZERO_PATTERN,&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
1237       }
1238       if (!sub_schurs->sum_S_Ej_inv_all && deluxe) {
1239         ierr = MatDuplicate(sub_schurs->S_Ej_all,MAT_SHARE_NONZERO_PATTERN,&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
1240       }
1241     }
1242 
1243     /* S_Ej_all */
1244     cum = cum2 = 0;
1245     ierr = MatDenseGetArray(S_all,&S_data);CHKERRQ(ierr);
1246     for (i=0;i<sub_schurs->n_subs;i++) {
1247       PetscInt j;
1248 
1249       /* get S_E */
1250       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
1251       if (S_lower_triangular) { /* I need to expand the upper triangular data (column oriented) */
1252         PetscInt k;
1253         for (k=0;k<subset_size;k++) {
1254           for (j=k;j<subset_size;j++) {
1255             work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
1256             work[j*subset_size+k] = PetscConj(S_data[cum2+k*size_schur+j]);
1257           }
1258         }
1259       } else { /* just copy to workspace */
1260         PetscInt k;
1261         for (k=0;k<subset_size;k++) {
1262           for (j=0;j<subset_size;j++) {
1263             work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
1264           }
1265         }
1266       }
1267       /* insert S_E values */
1268       for (j=0;j<subset_size;j++) dummy_idx[j] = cum+j;
1269       if (sub_schurs->change) {
1270         Mat change_sub,SEj,T;
1271 
1272         /* change basis */
1273         ierr = KSPGetOperators(sub_schurs->change[i],&change_sub,NULL);CHKERRQ(ierr);
1274         ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&SEj);CHKERRQ(ierr);
1275         if (!sub_schurs->change_with_qr) { /* currently there's no support for PtAP with P SeqAIJ */
1276           Mat T2;
1277           ierr = MatTransposeMatMult(change_sub,SEj,MAT_INITIAL_MATRIX,1.0,&T2);CHKERRQ(ierr);
1278           ierr = MatMatMult(T2,change_sub,MAT_INITIAL_MATRIX,1.0,&T);CHKERRQ(ierr);
1279           ierr = MatConvert(T,MATSEQDENSE,MAT_INPLACE_MATRIX,&T);CHKERRQ(ierr);
1280           ierr = MatDestroy(&T2);CHKERRQ(ierr);
1281         } else {
1282           ierr = MatPtAP(SEj,change_sub,MAT_INITIAL_MATRIX,1.0,&T);CHKERRQ(ierr);
1283         }
1284         ierr = MatCopy(T,SEj,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
1285         ierr = MatDestroy(&T);CHKERRQ(ierr);
1286         ierr = MatZeroRowsColumnsIS(SEj,sub_schurs->change_primal_sub[i],1.0,NULL,NULL);CHKERRQ(ierr);
1287         ierr = MatDestroy(&SEj);CHKERRQ(ierr);
1288       }
1289       if (deluxe) {
1290         ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
1291         /* if adaptivity is requested, invert S_E blocks */
1292         if (compute_Stilda) {
1293           PetscInt k;
1294 
1295           ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr);
1296           ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
1297           if (use_potr) {
1298             PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,work,&B_N,&B_ierr));
1299             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
1300             PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,work,&B_N,&B_ierr));
1301             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
1302             for (k=0;k<subset_size;k++) {
1303               for (j=k;j<subset_size;j++) {
1304                 work[j*subset_size+k] = work[k*subset_size+j];
1305               }
1306             }
1307           } else if (use_sytr) {
1308             PetscStackCallBLAS("LAPACKsytrf",LAPACKsytrf_("L",&B_N,work,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
1309             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYTRF Lapack routine %d",(int)B_ierr);
1310             PetscStackCallBLAS("LAPACKsytri",LAPACKsytri_("L",&B_N,work,&B_N,pivots,Bwork,&B_ierr));
1311             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYTRI Lapack routine %d",(int)B_ierr);
1312             for (k=0;k<subset_size;k++) {
1313               for (j=k;j<subset_size;j++) {
1314                 work[j*subset_size+k] = work[k*subset_size+j];
1315               }
1316             }
1317           } else {
1318             PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,work,&B_N,pivots,&B_ierr));
1319             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
1320             PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,work,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
1321             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
1322           }
1323           ierr = PetscLogFlops(1.0*subset_size*subset_size*subset_size);CHKERRQ(ierr);
1324           ierr = PetscFPTrapPop();CHKERRQ(ierr);
1325           ierr = MatSetValues(sub_schurs->sum_S_Ej_inv_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
1326         }
1327       } else if (compute_Stilda) { /* not using deluxe */
1328         Mat         SEj;
1329         Vec         D;
1330         PetscScalar *array;
1331 
1332         ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&SEj);CHKERRQ(ierr);
1333         ierr = VecGetArray(Dall,&array);CHKERRQ(ierr);
1334         ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,subset_size,array+cum,&D);CHKERRQ(ierr);
1335         ierr = VecRestoreArray(Dall,&array);CHKERRQ(ierr);
1336         ierr = VecShift(D,-1.);CHKERRQ(ierr);
1337         ierr = MatDiagonalScale(SEj,D,D);CHKERRQ(ierr);
1338         ierr = MatDestroy(&SEj);CHKERRQ(ierr);
1339         ierr = VecDestroy(&D);CHKERRQ(ierr);
1340         ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
1341       }
1342       cum += subset_size;
1343       cum2 += subset_size*(size_schur + 1);
1344     }
1345     ierr = MatDenseRestoreArray(S_all,&S_data);CHKERRQ(ierr);
1346 
1347     if (solver_S) {
1348       ierr = MatFactorRestoreSchurComplement(F,&S_all,MAT_FACTOR_SCHUR_UNFACTORED);CHKERRQ(ierr);
1349     }
1350 
1351     schur_factor = NULL;
1352     if (compute_Stilda && size_active_schur) {
1353 
1354       if (sub_schurs->n_subs == 1 && size_schur == size_active_schur && deluxe) { /* we already computed the inverse */
1355         PetscInt j;
1356         for (j=0;j<size_schur;j++) dummy_idx[j] = j;
1357         ierr = MatSetValues(sub_schurs->sum_S_Ej_tilda_all,size_schur,dummy_idx,size_schur,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
1358       } else {
1359         Mat S_all_inv=NULL;
1360         if (solver_S) {
1361           /* for adaptive selection we need S^-1; for solver reusage we need S_\Delta\Delta^-1.
1362              The latter is not the principal subminor for S^-1. However, the factors can be reused since S_\Delta\Delta is the leading principal submatrix of S */
1363           if (factor_workaround) {/* invert without calling MatFactorInvertSchurComplement, since we are hacking */
1364             PetscScalar *data;
1365             PetscInt     nd = 0;
1366 
1367             if (!use_potr) {
1368               SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor update not yet implemented for non SPD matrices");
1369             }
1370             ierr = MatFactorGetSchurComplement(F,&S_all_inv,NULL);CHKERRQ(ierr);
1371             ierr = MatDenseGetArray(S_all_inv,&data);CHKERRQ(ierr);
1372             if (sub_schurs->is_dir) { /* dirichlet dofs could have different scalings */
1373               ierr = ISGetLocalSize(sub_schurs->is_dir,&nd);CHKERRQ(ierr);
1374             }
1375 
1376             /* factor and invert activedofs and vertices (dirichlet dofs does not contribute) */
1377             if (schur_has_vertices) {
1378               Mat          M;
1379               PetscScalar *tdata;
1380               PetscInt     nv = 0, news;
1381 
1382               ierr = ISGetLocalSize(sub_schurs->is_vertices,&nv);CHKERRQ(ierr);
1383               news = size_active_schur + nv;
1384               ierr = PetscCalloc1(news*news,&tdata);CHKERRQ(ierr);
1385               for (i=0;i<size_active_schur;i++) {
1386                 ierr = PetscMemcpy(tdata+i*(news+1),data+i*(size_schur+1),(size_active_schur-i)*sizeof(PetscScalar));CHKERRQ(ierr);
1387                 ierr = PetscMemcpy(tdata+i*(news+1)+size_active_schur-i,data+i*size_schur+size_active_schur+nd,nv*sizeof(PetscScalar));CHKERRQ(ierr);
1388               }
1389               for (i=0;i<nv;i++) {
1390                 PetscInt k = i+size_active_schur;
1391                 ierr = PetscMemcpy(tdata+k*(news+1),data+(k+nd)*(size_schur+1),(nv-i)*sizeof(PetscScalar));CHKERRQ(ierr);
1392               }
1393 
1394               ierr = MatCreateSeqDense(PETSC_COMM_SELF,news,news,tdata,&M);CHKERRQ(ierr);
1395               ierr = MatSetOption(M,MAT_SPD,PETSC_TRUE);CHKERRQ(ierr);
1396               ierr = MatCholeskyFactor(M,NULL,NULL);CHKERRQ(ierr);
1397               /* save the factors */
1398               cum = 0;
1399               ierr = PetscMalloc1((size_active_schur*(size_active_schur +1))/2+nd,&schur_factor);CHKERRQ(ierr);
1400               for (i=0;i<size_active_schur;i++) {
1401                 ierr = PetscMemcpy(schur_factor+cum,tdata+i*(news+1),(size_active_schur-i)*sizeof(PetscScalar));CHKERRQ(ierr);
1402                 cum += size_active_schur - i;
1403               }
1404               for (i=0;i<nd;i++) schur_factor[cum+i] = PetscSqrtReal(PetscRealPart(data[(i+size_active_schur)*(size_schur+1)]));
1405               ierr = MatSeqDenseInvertFactors_Private(M);CHKERRQ(ierr);
1406               /* move back just the active dofs to the Schur complement */
1407               for (i=0;i<size_active_schur;i++) {
1408                 ierr = PetscMemcpy(data+i*size_schur,tdata+i*news,size_active_schur*sizeof(PetscScalar));CHKERRQ(ierr);
1409               }
1410               ierr = PetscFree(tdata);CHKERRQ(ierr);
1411               ierr = MatDestroy(&M);CHKERRQ(ierr);
1412             } else { /* we can factorize and invert just the activedofs */
1413               Mat         M;
1414               PetscScalar *aux;
1415 
1416               ierr = PetscMalloc1(nd,&aux);CHKERRQ(ierr);
1417               for (i=0;i<nd;i++) aux[i] = 1.0/data[(i+size_active_schur)*(size_schur+1)];
1418               ierr = MatCreateSeqDense(PETSC_COMM_SELF,size_active_schur,size_active_schur,data,&M);CHKERRQ(ierr);
1419               ierr = MatSeqDenseSetLDA(M,size_schur);CHKERRQ(ierr);
1420               ierr = MatSetOption(M,MAT_SPD,PETSC_TRUE);CHKERRQ(ierr);
1421               ierr = MatCholeskyFactor(M,NULL,NULL);CHKERRQ(ierr);
1422               ierr = MatSeqDenseInvertFactors_Private(M);CHKERRQ(ierr);
1423               ierr = MatDestroy(&M);CHKERRQ(ierr);
1424               ierr = MatCreateSeqDense(PETSC_COMM_SELF,size_schur,nd,data+size_active_schur*size_schur,&M);CHKERRQ(ierr);
1425               ierr = MatZeroEntries(M);CHKERRQ(ierr);
1426               ierr = MatDestroy(&M);CHKERRQ(ierr);
1427               ierr = MatCreateSeqDense(PETSC_COMM_SELF,nd,size_schur,data+size_active_schur,&M);CHKERRQ(ierr);
1428               ierr = MatSeqDenseSetLDA(M,size_schur);CHKERRQ(ierr);
1429               ierr = MatZeroEntries(M);CHKERRQ(ierr);
1430               ierr = MatDestroy(&M);CHKERRQ(ierr);
1431               for (i=0;i<nd;i++) data[(i+size_active_schur)*(size_schur+1)] = aux[i];
1432               ierr = PetscFree(aux);CHKERRQ(ierr);
1433             }
1434             ierr = MatDenseRestoreArray(S_all_inv,&data);CHKERRQ(ierr);
1435           } else { /* use MatFactor calls to invert S */
1436             ierr = MatFactorInvertSchurComplement(F);CHKERRQ(ierr);
1437             ierr = MatFactorGetSchurComplement(F,&S_all_inv,NULL);CHKERRQ(ierr);
1438           }
1439         } else { /* we need to invert explicitly since we are not using MatFactor for S */
1440           ierr = PetscObjectReference((PetscObject)S_all);CHKERRQ(ierr);
1441           S_all_inv = S_all;
1442           ierr = MatDenseGetArray(S_all_inv,&S_data);CHKERRQ(ierr);
1443           ierr = PetscBLASIntCast(size_schur,&B_N);CHKERRQ(ierr);
1444           ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
1445           if (use_potr) {
1446             PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,S_data,&B_N,&B_ierr));
1447             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
1448             PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,S_data,&B_N,&B_ierr));
1449             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
1450           } else if (use_sytr) {
1451             PetscStackCallBLAS("LAPACKsytrf",LAPACKsytrf_("L",&B_N,S_data,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
1452             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYTRF Lapack routine %d",(int)B_ierr);
1453             PetscStackCallBLAS("LAPACKsytri",LAPACKsytri_("L",&B_N,S_data,&B_N,pivots,Bwork,&B_ierr));
1454             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYTRI Lapack routine %d",(int)B_ierr);
1455           } else {
1456             PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,S_data,&B_N,pivots,&B_ierr));
1457             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
1458             PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,S_data,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
1459             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
1460           }
1461           ierr = PetscLogFlops(1.0*size_schur*size_schur*size_schur);CHKERRQ(ierr);
1462           ierr = PetscFPTrapPop();CHKERRQ(ierr);
1463           ierr = MatDenseRestoreArray(S_all_inv,&S_data);CHKERRQ(ierr);
1464         }
1465         /* S_Ej_tilda_all */
1466         cum = cum2 = 0;
1467         ierr = MatDenseGetArray(S_all_inv,&S_data);CHKERRQ(ierr);
1468         for (i=0;i<sub_schurs->n_subs;i++) {
1469           PetscInt j;
1470 
1471           ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
1472           /* get (St^-1)_E */
1473           /* Unless we are changing the variables, I don't need to expand to upper triangular since St^-1
1474              will be properly accessed later during adaptive selection */
1475           if (S_lower_triangular) {
1476             PetscInt k;
1477             if (sub_schurs->change) {
1478               for (k=0;k<subset_size;k++) {
1479                 for (j=k;j<subset_size;j++) {
1480                   work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
1481                   work[j*subset_size+k] = work[k*subset_size+j];
1482                 }
1483               }
1484             } else {
1485               for (k=0;k<subset_size;k++) {
1486                 for (j=k;j<subset_size;j++) {
1487                   work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
1488                 }
1489               }
1490             }
1491           } else {
1492             PetscInt k;
1493             for (k=0;k<subset_size;k++) {
1494               for (j=0;j<subset_size;j++) {
1495                 work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
1496               }
1497             }
1498           }
1499           if (sub_schurs->change) {
1500             Mat change_sub,SEj,T;
1501 
1502             /* change basis */
1503             ierr = KSPGetOperators(sub_schurs->change[i],&change_sub,NULL);CHKERRQ(ierr);
1504             ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&SEj);CHKERRQ(ierr);
1505             if (!sub_schurs->change_with_qr) { /* currently there's no support for PtAP with P SeqAIJ */
1506               Mat T2;
1507               ierr = MatTransposeMatMult(change_sub,SEj,MAT_INITIAL_MATRIX,1.0,&T2);CHKERRQ(ierr);
1508               ierr = MatMatMult(T2,change_sub,MAT_INITIAL_MATRIX,1.0,&T);CHKERRQ(ierr);
1509               ierr = MatDestroy(&T2);CHKERRQ(ierr);
1510               ierr = MatConvert(T,MATSEQDENSE,MAT_INPLACE_MATRIX,&T);CHKERRQ(ierr);
1511             } else {
1512               ierr = MatPtAP(SEj,change_sub,MAT_INITIAL_MATRIX,1.0,&T);CHKERRQ(ierr);
1513             }
1514             ierr = MatCopy(T,SEj,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
1515             ierr = MatDestroy(&T);CHKERRQ(ierr);
1516             /* set diagonal entry to a very large value to pick the basis we are eliminating as the first eigenvectors with adaptive selection */
1517             ierr = MatZeroRowsColumnsIS(SEj,sub_schurs->change_primal_sub[i],1./PETSC_SMALL,NULL,NULL);CHKERRQ(ierr);
1518             ierr = MatDestroy(&SEj);CHKERRQ(ierr);
1519           }
1520           for (j=0;j<subset_size;j++) dummy_idx[j] = cum+j;
1521           ierr = MatSetValues(sub_schurs->sum_S_Ej_tilda_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
1522           cum += subset_size;
1523           cum2 += subset_size*(size_schur + 1);
1524         }
1525         ierr = MatDenseRestoreArray(S_all_inv,&S_data);CHKERRQ(ierr);
1526         if (solver_S) {
1527           if (schur_has_vertices) {
1528             ierr = MatFactorRestoreSchurComplement(F,&S_all_inv,MAT_FACTOR_SCHUR_FACTORED);CHKERRQ(ierr);
1529           } else {
1530             ierr = MatFactorRestoreSchurComplement(F,&S_all_inv,MAT_FACTOR_SCHUR_INVERTED);CHKERRQ(ierr);
1531           }
1532         }
1533         ierr = MatDestroy(&S_all_inv);CHKERRQ(ierr);
1534       }
1535 
1536       /* move back factors if needed */
1537       if (schur_has_vertices) {
1538         Mat      S_tmp;
1539         PetscInt nd = 0;
1540 
1541         if (!solver_S) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"This should not happen");
1542         ierr = MatFactorGetSchurComplement(F,&S_tmp,NULL);CHKERRQ(ierr);
1543         if (use_potr) {
1544           PetscScalar *data;
1545 
1546           ierr = MatDenseGetArray(S_tmp,&data);CHKERRQ(ierr);
1547           ierr = PetscMemzero(data,size_schur*size_schur*sizeof(PetscScalar));CHKERRQ(ierr);
1548 
1549           if (S_lower_triangular) {
1550             cum = 0;
1551             for (i=0;i<size_active_schur;i++) {
1552               ierr = PetscMemcpy(data+i*(size_schur+1),schur_factor+cum,(size_active_schur-i)*sizeof(PetscScalar));CHKERRQ(ierr);
1553               cum += size_active_schur-i;
1554 	    }
1555           } else {
1556             ierr = PetscMemcpy(data,schur_factor,size_schur*size_schur*sizeof(PetscScalar));CHKERRQ(ierr);
1557           }
1558           if (sub_schurs->is_dir) {
1559             ierr = ISGetLocalSize(sub_schurs->is_dir,&nd);CHKERRQ(ierr);
1560 	    for (i=0;i<nd;i++) {
1561 	      data[(i+size_active_schur)*(size_schur+1)] = schur_factor[cum+i];
1562 	    }
1563           }
1564           /* workaround: since I cannot modify the matrices used inside the solvers for the forward and backward substitutions,
1565              set the diagonal entry of the Schur factor to a very large value */
1566           for (i=size_active_schur+nd;i<size_schur;i++) {
1567             data[i*(size_schur+1)] = infty;
1568           }
1569           ierr = MatDenseRestoreArray(S_tmp,&data);CHKERRQ(ierr);
1570         } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor update not yet implemented for non SPD matrices");
1571         ierr = MatFactorRestoreSchurComplement(F,&S_tmp,MAT_FACTOR_SCHUR_FACTORED);CHKERRQ(ierr);
1572       }
1573     } else if (factor_workaround) { /* we need to eliminate any unneeded coupling */
1574       PetscScalar *data;
1575       PetscInt    nd = 0;
1576 
1577       if (sub_schurs->is_dir) { /* dirichlet dofs could have different scalings */
1578         ierr = ISGetLocalSize(sub_schurs->is_dir,&nd);CHKERRQ(ierr);
1579       }
1580       ierr = MatFactorGetSchurComplement(F,&S_all,NULL);CHKERRQ(ierr);
1581       ierr = MatDenseGetArray(S_all,&data);CHKERRQ(ierr);
1582       for (i=0;i<size_active_schur;i++) {
1583         ierr = PetscMemzero(data+i*size_schur+size_active_schur,(size_schur-size_active_schur)*sizeof(PetscScalar));CHKERRQ(ierr);
1584       }
1585       for (i=size_active_schur+nd;i<size_schur;i++) {
1586         ierr = PetscMemzero(data+i*size_schur+size_active_schur,(size_schur-size_active_schur)*sizeof(PetscScalar));CHKERRQ(ierr);
1587         data[i*(size_schur+1)] = infty;
1588       }
1589       ierr = MatDenseRestoreArray(S_all,&data);CHKERRQ(ierr);
1590       ierr = MatFactorRestoreSchurComplement(F,&S_all,MAT_FACTOR_SCHUR_UNFACTORED);CHKERRQ(ierr);
1591     }
1592     ierr = PetscFree2(dummy_idx,work);CHKERRQ(ierr);
1593     ierr = PetscFree(schur_factor);CHKERRQ(ierr);
1594     ierr = VecDestroy(&Dall);CHKERRQ(ierr);
1595   }
1596   ierr = PetscFree(nnz);CHKERRQ(ierr);
1597   ierr = MatDestroy(&F);CHKERRQ(ierr);
1598   ierr = ISDestroy(&is_I_layer);CHKERRQ(ierr);
1599   ierr = MatDestroy(&S_all);CHKERRQ(ierr);
1600   ierr = MatDestroy(&A_BB);CHKERRQ(ierr);
1601   ierr = MatDestroy(&A_IB);CHKERRQ(ierr);
1602   ierr = MatDestroy(&A_BI);CHKERRQ(ierr);
1603 
1604   /* speed up matrix assembly */
1605   ierr = PetscMalloc1(sub_schurs->n_subs,&nnz);CHKERRQ(ierr);
1606   for (i=0;i<sub_schurs->n_subs;i++) {
1607     ierr = ISGetLocalSize(sub_schurs->is_subs[i],&nnz[i]);CHKERRQ(ierr);
1608   }
1609   ierr = ISCreateGeneral(PETSC_COMM_SELF,sub_schurs->n_subs,nnz,PETSC_OWN_POINTER,&is_I_layer);CHKERRQ(ierr);
1610   ierr = MatSetVariableBlockSizes(sub_schurs->S_Ej_all,sub_schurs->n_subs,nnz);CHKERRQ(ierr);
1611   ierr = MatAssemblyBegin(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1612   ierr = MatAssemblyEnd(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1613   if (compute_Stilda) {
1614     ierr = MatSetVariableBlockSizes(sub_schurs->sum_S_Ej_tilda_all,sub_schurs->n_subs,nnz);CHKERRQ(ierr);
1615     ierr = MatAssemblyBegin(sub_schurs->sum_S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1616     ierr = MatAssemblyEnd(sub_schurs->sum_S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1617     if (deluxe) {
1618       ierr = MatSetVariableBlockSizes(sub_schurs->sum_S_Ej_inv_all,sub_schurs->n_subs,nnz);CHKERRQ(ierr);
1619       ierr = MatAssemblyBegin(sub_schurs->sum_S_Ej_inv_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1620       ierr = MatAssemblyEnd(sub_schurs->sum_S_Ej_inv_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1621     }
1622   }
1623   ierr = ISDestroy(&is_I_layer);CHKERRQ(ierr);
1624 
1625   /* Global matrix of all assembled Schur on subsets */
1626   ierr = MatISSetLocalMat(work_mat,sub_schurs->S_Ej_all);CHKERRQ(ierr);
1627   ierr = MatAssemblyBegin(work_mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1628   ierr = MatAssemblyEnd(work_mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1629   ierr = MatISSetMPIXAIJPreallocation_Private(work_mat,global_schur_subsets,PETSC_TRUE);CHKERRQ(ierr);
1630   ierr = MatConvert(work_mat,MATAIJ,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
1631 
1632   /* Get local part of (\sum_j S_Ej) */
1633   ierr = MatCreateSubMatrices(global_schur_subsets,1,&all_subsets_n,&all_subsets_n,MAT_INITIAL_MATRIX,&submats);CHKERRQ(ierr);
1634   if (!sub_schurs->sum_S_Ej_all) {
1635     ierr = MatDuplicate(submats[0],MAT_COPY_VALUES,&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
1636   } else {
1637     ierr = MatCopy(submats[0],sub_schurs->sum_S_Ej_all,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
1638   }
1639 
1640   /* Get local part of (\sum_j S^-1_Ej) (\sum_j St^-1_Ej) */
1641   if (compute_Stilda) {
1642     ierr = MatISSetLocalMat(work_mat,sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
1643     ierr = MatConvert(work_mat,MATAIJ,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
1644     ierr = MatCreateSubMatrices(global_schur_subsets,1,&all_subsets_n,&all_subsets_n,MAT_REUSE_MATRIX,&submats);CHKERRQ(ierr);
1645     ierr = MatCopy(submats[0],sub_schurs->sum_S_Ej_tilda_all,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
1646     if (deluxe) {
1647       ierr = MatISSetLocalMat(work_mat,sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
1648       ierr = MatConvert(work_mat,MATAIJ,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
1649       ierr = MatCreateSubMatrices(global_schur_subsets,1,&all_subsets_n,&all_subsets_n,MAT_REUSE_MATRIX,&submats);CHKERRQ(ierr);
1650       ierr = MatCopy(submats[0],sub_schurs->sum_S_Ej_inv_all,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
1651     } else {
1652       PetscScalar *array;
1653       PetscInt    cum;
1654 
1655       ierr = MatSeqAIJGetArray(sub_schurs->sum_S_Ej_tilda_all,&array);CHKERRQ(ierr);
1656       cum = 0;
1657       for (i=0;i<sub_schurs->n_subs;i++) {
1658         ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
1659         ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr);
1660         ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
1661         if (use_potr) {
1662           PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,array+cum,&B_N,&B_ierr));
1663           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
1664           PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,array+cum,&B_N,&B_ierr));
1665           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
1666         } else if (use_sytr) {
1667           PetscStackCallBLAS("LAPACKsytrf",LAPACKsytrf_("L",&B_N,array+cum,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
1668           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYTRF Lapack routine %d",(int)B_ierr);
1669           PetscStackCallBLAS("LAPACKsytri",LAPACKsytri_("L",&B_N,array+cum,&B_N,pivots,Bwork,&B_ierr));
1670           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYTRI Lapack routine %d",(int)B_ierr);
1671         } else {
1672           PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,array+cum,&B_N,pivots,&B_ierr));
1673           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
1674           PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,array+cum,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
1675           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
1676         }
1677         ierr = PetscLogFlops(1.0*subset_size*subset_size*subset_size);CHKERRQ(ierr);
1678         ierr = PetscFPTrapPop();CHKERRQ(ierr);
1679         cum += subset_size*subset_size;
1680       }
1681       ierr = MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_tilda_all,&array);CHKERRQ(ierr);
1682       ierr = PetscObjectReference((PetscObject)sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
1683       ierr = MatDestroy(&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
1684       sub_schurs->sum_S_Ej_inv_all = sub_schurs->sum_S_Ej_all;
1685     }
1686   }
1687   ierr = MatDestroySubMatrices(1,&submats);CHKERRQ(ierr);
1688   if (print_schurs) {
1689     PetscViewer viewer;
1690     char        filename[256];
1691     PetscInt    cum;
1692 
1693     ierr = PetscSNPrintf(filename,sizeof(filename),"sub_schurs_mats_r%d.m",PetscGlobalRank);CHKERRQ(ierr);
1694     ierr = PetscViewerASCIIOpen(PETSC_COMM_SELF,filename,&viewer);CHKERRQ(ierr);
1695     ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);CHKERRQ(ierr);
1696     if (sub_schurs->S_Ej_all) {
1697       ierr = PetscObjectSetName((PetscObject)sub_schurs->S_Ej_all,"SE");CHKERRQ(ierr);
1698       ierr = MatView(sub_schurs->S_Ej_all,viewer);CHKERRQ(ierr);
1699     }
1700     if (sub_schurs->sum_S_Ej_all) {
1701       ierr = PetscObjectSetName((PetscObject)sub_schurs->sum_S_Ej_all,"SSE");CHKERRQ(ierr);
1702       ierr = MatView(sub_schurs->sum_S_Ej_all,viewer);CHKERRQ(ierr);
1703     }
1704     if (sub_schurs->sum_S_Ej_inv_all) {
1705       ierr = PetscObjectSetName((PetscObject)sub_schurs->sum_S_Ej_inv_all,"SSEm");CHKERRQ(ierr);
1706       ierr = MatView(sub_schurs->sum_S_Ej_inv_all,viewer);CHKERRQ(ierr);
1707     }
1708     if (sub_schurs->sum_S_Ej_tilda_all) {
1709       ierr = PetscObjectSetName((PetscObject)sub_schurs->sum_S_Ej_tilda_all,"SSEt");CHKERRQ(ierr);
1710       ierr = MatView(sub_schurs->sum_S_Ej_tilda_all,viewer);CHKERRQ(ierr);
1711     }
1712     for (i=0,cum=0;i<sub_schurs->n_subs;i++) {
1713       IS   is;
1714       char name[16];
1715 
1716       ierr = PetscSNPrintf(name,sizeof(name),"IE%D",i);CHKERRQ(ierr);
1717       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
1718       ierr = ISCreateStride(PETSC_COMM_SELF,subset_size,cum,1,&is);CHKERRQ(ierr);
1719       ierr = PetscObjectSetName((PetscObject)is,name);CHKERRQ(ierr);
1720       ierr = ISView(is,viewer);CHKERRQ(ierr);
1721       ierr = ISDestroy(&is);CHKERRQ(ierr);
1722       cum += subset_size;
1723     }
1724     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
1725   }
1726 
1727   /* free workspace */
1728   ierr = PetscFree(submats);CHKERRQ(ierr);
1729   ierr = PetscFree2(Bwork,pivots);CHKERRQ(ierr);
1730   ierr = MatDestroy(&global_schur_subsets);CHKERRQ(ierr);
1731   ierr = MatDestroy(&work_mat);CHKERRQ(ierr);
1732   ierr = ISDestroy(&all_subsets_n);CHKERRQ(ierr);
1733   ierr = PetscCommDestroy(&comm_n);CHKERRQ(ierr);
1734   PetscFunctionReturn(0);
1735 }
1736 
1737 PetscErrorCode PCBDDCSubSchursInit(PCBDDCSubSchurs sub_schurs, const char* prefix, IS is_I, IS is_B, PCBDDCGraph graph, ISLocalToGlobalMapping BtoNmap, PetscBool copycc)
1738 {
1739   IS              *faces,*edges,*all_cc,vertices;
1740   PetscInt        i,n_faces,n_edges,n_all_cc;
1741   PetscBool       is_sorted,ispetsc;
1742   PetscErrorCode  ierr;
1743 
1744   PetscFunctionBegin;
1745   ierr = ISSorted(is_I,&is_sorted);CHKERRQ(ierr);
1746   if (!is_sorted) SETERRQ(PetscObjectComm((PetscObject)is_I),PETSC_ERR_PLIB,"IS for I dofs should be shorted");
1747   ierr = ISSorted(is_B,&is_sorted);CHKERRQ(ierr);
1748   if (!is_sorted) SETERRQ(PetscObjectComm((PetscObject)is_B),PETSC_ERR_PLIB,"IS for B dofs should be shorted");
1749 
1750   /* reset any previous data */
1751   ierr = PCBDDCSubSchursReset(sub_schurs);CHKERRQ(ierr);
1752 
1753   /* get index sets for faces and edges (already sorted by global ordering) */
1754   ierr = PCBDDCGraphGetCandidatesIS(graph,&n_faces,&faces,&n_edges,&edges,&vertices);CHKERRQ(ierr);
1755   n_all_cc = n_faces+n_edges;
1756   ierr = PetscBTCreate(n_all_cc,&sub_schurs->is_edge);CHKERRQ(ierr);
1757   ierr = PetscMalloc1(n_all_cc,&all_cc);CHKERRQ(ierr);
1758   for (i=0;i<n_faces;i++) {
1759     if (copycc) {
1760       ierr = ISDuplicate(faces[i],&all_cc[i]);CHKERRQ(ierr);
1761     } else {
1762       ierr = PetscObjectReference((PetscObject)faces[i]);CHKERRQ(ierr);
1763       all_cc[i] = faces[i];
1764     }
1765   }
1766   for (i=0;i<n_edges;i++) {
1767     if (copycc) {
1768       ierr = ISDuplicate(edges[i],&all_cc[n_faces+i]);CHKERRQ(ierr);
1769     } else {
1770       ierr = PetscObjectReference((PetscObject)edges[i]);CHKERRQ(ierr);
1771       all_cc[n_faces+i] = edges[i];
1772     }
1773     ierr = PetscBTSet(sub_schurs->is_edge,n_faces+i);CHKERRQ(ierr);
1774   }
1775   ierr = PetscObjectReference((PetscObject)vertices);CHKERRQ(ierr);
1776   sub_schurs->is_vertices = vertices;
1777   ierr = PCBDDCGraphRestoreCandidatesIS(graph,&n_faces,&faces,&n_edges,&edges,&vertices);CHKERRQ(ierr);
1778   sub_schurs->is_dir = NULL;
1779   ierr = PCBDDCGraphGetDirichletDofsB(graph,&sub_schurs->is_dir);CHKERRQ(ierr);
1780 
1781   /* Determine if MatFactor can be used */
1782   ierr = PetscStrallocpy(prefix,&sub_schurs->prefix);CHKERRQ(ierr);
1783 #if defined(PETSC_HAVE_MUMPS)
1784   ierr = PetscStrncpy(sub_schurs->mat_solver_type,MATSOLVERMUMPS,64);CHKERRQ(ierr);
1785 #elif defined(PETSC_HAVE_MKL_PARDISO)
1786   ierr = PetscStrncpy(sub_schurs->mat_solver_type,MATSOLVERMKL_PARDISO,64);CHKERRQ(ierr);
1787 #else
1788   ierr = PetscStrncpy(sub_schurs->mat_solver_type,MATSOLVERPETSC,64);CHKERRQ(ierr);
1789 #endif
1790 #if defined(PETSC_USE_COMPLEX)
1791   sub_schurs->is_hermitian  = PETSC_FALSE; /* Hermitian Cholesky is not supported by PETSc and external packages */
1792 #else
1793   sub_schurs->is_hermitian  = PETSC_TRUE;
1794 #endif
1795   sub_schurs->is_posdef     = PETSC_TRUE;
1796   sub_schurs->is_symmetric  = PETSC_TRUE;
1797   sub_schurs->debug         = PETSC_FALSE;
1798   sub_schurs->restrict_comm = PETSC_FALSE;
1799   ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)graph->l2gmap),sub_schurs->prefix,"BDDC sub_schurs options","PC");CHKERRQ(ierr);
1800   ierr = PetscOptionsString("-sub_schurs_mat_solver_type","Specific direct solver to use",NULL,sub_schurs->mat_solver_type,sub_schurs->mat_solver_type,64,NULL);CHKERRQ(ierr);
1801   ierr = PetscOptionsBool("-sub_schurs_symmetric","Symmetric problem",NULL,sub_schurs->is_symmetric,&sub_schurs->is_symmetric,NULL);CHKERRQ(ierr);
1802   ierr = PetscOptionsBool("-sub_schurs_hermitian","Hermitian problem",NULL,sub_schurs->is_hermitian,&sub_schurs->is_hermitian,NULL);CHKERRQ(ierr);
1803   ierr = PetscOptionsBool("-sub_schurs_posdef","Positive definite problem",NULL,sub_schurs->is_posdef,&sub_schurs->is_posdef,NULL);CHKERRQ(ierr);
1804   ierr = PetscOptionsBool("-sub_schurs_restrictcomm","Restrict communicator on active processes only",NULL,sub_schurs->restrict_comm,&sub_schurs->restrict_comm,NULL);CHKERRQ(ierr);
1805   ierr = PetscOptionsBool("-sub_schurs_debug","Debug output",NULL,sub_schurs->debug,&sub_schurs->debug,NULL);CHKERRQ(ierr);
1806   ierr = PetscOptionsEnd();CHKERRQ(ierr);
1807   ierr = PetscStrcmp(sub_schurs->mat_solver_type,MATSOLVERPETSC,&ispetsc);CHKERRQ(ierr);
1808   sub_schurs->schur_explicit = (PetscBool)!ispetsc;
1809 
1810   /* for reals, symmetric and hermitian are synonims */
1811 #if !defined(PETSC_USE_COMPLEX)
1812   sub_schurs->is_symmetric = (PetscBool)(sub_schurs->is_symmetric && sub_schurs->is_hermitian);
1813   sub_schurs->is_hermitian = sub_schurs->is_symmetric;
1814 #endif
1815 
1816   ierr = PetscObjectReference((PetscObject)is_I);CHKERRQ(ierr);
1817   sub_schurs->is_I = is_I;
1818   ierr = PetscObjectReference((PetscObject)is_B);CHKERRQ(ierr);
1819   sub_schurs->is_B = is_B;
1820   ierr = PetscObjectReference((PetscObject)graph->l2gmap);CHKERRQ(ierr);
1821   sub_schurs->l2gmap = graph->l2gmap;
1822   ierr = PetscObjectReference((PetscObject)BtoNmap);CHKERRQ(ierr);
1823   sub_schurs->BtoNmap = BtoNmap;
1824   sub_schurs->n_subs = n_all_cc;
1825   sub_schurs->is_subs = all_cc;
1826   sub_schurs->S_Ej_all = NULL;
1827   sub_schurs->sum_S_Ej_all = NULL;
1828   sub_schurs->sum_S_Ej_inv_all = NULL;
1829   sub_schurs->sum_S_Ej_tilda_all = NULL;
1830   sub_schurs->is_Ej_all = NULL;
1831   PetscFunctionReturn(0);
1832 }
1833 
1834 PetscErrorCode PCBDDCSubSchursCreate(PCBDDCSubSchurs *sub_schurs)
1835 {
1836   PCBDDCSubSchurs schurs_ctx;
1837   PetscErrorCode  ierr;
1838 
1839   PetscFunctionBegin;
1840   ierr = PetscNew(&schurs_ctx);CHKERRQ(ierr);
1841   schurs_ctx->n_subs = 0;
1842   *sub_schurs = schurs_ctx;
1843   PetscFunctionReturn(0);
1844 }
1845 
1846 PetscErrorCode PCBDDCSubSchursReset(PCBDDCSubSchurs sub_schurs)
1847 {
1848   PetscInt       i;
1849   PetscErrorCode ierr;
1850 
1851   PetscFunctionBegin;
1852   if (!sub_schurs) PetscFunctionReturn(0);
1853   ierr = PetscFree(sub_schurs->prefix);CHKERRQ(ierr);
1854   ierr = MatDestroy(&sub_schurs->A);CHKERRQ(ierr);
1855   ierr = MatDestroy(&sub_schurs->S);CHKERRQ(ierr);
1856   ierr = ISDestroy(&sub_schurs->is_I);CHKERRQ(ierr);
1857   ierr = ISDestroy(&sub_schurs->is_B);CHKERRQ(ierr);
1858   ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->l2gmap);CHKERRQ(ierr);
1859   ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->BtoNmap);CHKERRQ(ierr);
1860   ierr = MatDestroy(&sub_schurs->S_Ej_all);CHKERRQ(ierr);
1861   ierr = MatDestroy(&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
1862   ierr = MatDestroy(&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
1863   ierr = MatDestroy(&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
1864   ierr = ISDestroy(&sub_schurs->is_Ej_all);CHKERRQ(ierr);
1865   ierr = ISDestroy(&sub_schurs->is_vertices);CHKERRQ(ierr);
1866   ierr = ISDestroy(&sub_schurs->is_dir);CHKERRQ(ierr);
1867   ierr = PetscBTDestroy(&sub_schurs->is_edge);CHKERRQ(ierr);
1868   for (i=0;i<sub_schurs->n_subs;i++) {
1869     ierr = ISDestroy(&sub_schurs->is_subs[i]);CHKERRQ(ierr);
1870   }
1871   if (sub_schurs->n_subs) {
1872     ierr = PetscFree(sub_schurs->is_subs);CHKERRQ(ierr);
1873   }
1874   if (sub_schurs->reuse_solver) {
1875     ierr = PCBDDCReuseSolversReset(sub_schurs->reuse_solver);CHKERRQ(ierr);
1876   }
1877   ierr = PetscFree(sub_schurs->reuse_solver);CHKERRQ(ierr);
1878   if (sub_schurs->change) {
1879     for (i=0;i<sub_schurs->n_subs;i++) {
1880       ierr = KSPDestroy(&sub_schurs->change[i]);CHKERRQ(ierr);
1881       ierr = ISDestroy(&sub_schurs->change_primal_sub[i]);CHKERRQ(ierr);
1882     }
1883   }
1884   ierr = PetscFree(sub_schurs->change);CHKERRQ(ierr);
1885   ierr = PetscFree(sub_schurs->change_primal_sub);CHKERRQ(ierr);
1886   sub_schurs->n_subs = 0;
1887   PetscFunctionReturn(0);
1888 }
1889 
1890 PetscErrorCode PCBDDCSubSchursDestroy(PCBDDCSubSchurs* sub_schurs)
1891 {
1892   PetscErrorCode ierr;
1893 
1894   PetscFunctionBegin;
1895   ierr = PCBDDCSubSchursReset(*sub_schurs);CHKERRQ(ierr);
1896   ierr = PetscFree(*sub_schurs);CHKERRQ(ierr);
1897   PetscFunctionReturn(0);
1898 }
1899 
1900 PETSC_STATIC_INLINE PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt* queue_tip,PetscInt n_prev,PetscBT touched,PetscInt* xadj,PetscInt* adjncy,PetscInt* n_added)
1901 {
1902   PetscInt       i,j,n;
1903   PetscErrorCode ierr;
1904 
1905   PetscFunctionBegin;
1906   n = 0;
1907   for (i=-n_prev;i<0;i++) {
1908     PetscInt start_dof = queue_tip[i];
1909     for (j=xadj[start_dof];j<xadj[start_dof+1];j++) {
1910       PetscInt dof = adjncy[j];
1911       if (!PetscBTLookup(touched,dof)) {
1912         ierr = PetscBTSet(touched,dof);CHKERRQ(ierr);
1913         queue_tip[n] = dof;
1914         n++;
1915       }
1916     }
1917   }
1918   *n_added = n;
1919   PetscFunctionReturn(0);
1920 }
1921