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