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