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