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