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