xref: /petsc/src/ksp/pc/impls/bddc/bddcschurs.c (revision 91af69089d34f7e08ffee5b1dcb09643a4c10217)
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) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not for parallel matrices");
298   if (reuse == MAT_REUSE_MATRIX) {
299     PetscBool Sdense;
300 
301     ierr = PetscObjectTypeCompare((PetscObject)*S, MATSEQDENSE, &Sdense);CHKERRQ(ierr);
302     if (!Sdense) SETERRQ(PetscObjectComm((PetscObject)M),PETSC_ERR_SUP,"S should dense");
303   }
304   ierr = MatSchurComplementGetSubMatrices(M, NULL, NULL, &B, &C, &D);CHKERRQ(ierr);
305   ierr = MatSchurComplementGetKSP(M, &ksp);CHKERRQ(ierr);
306   ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr);
307   ierr = PetscObjectTypeCompare((PetscObject) pc, PCLU, &isLU);CHKERRQ(ierr);
308   ierr = PetscObjectTypeCompare((PetscObject) pc, PCILU, &isILU);CHKERRQ(ierr);
309   ierr = PetscObjectTypeCompare((PetscObject) pc, PCCHOLESKY, &isCHOL);CHKERRQ(ierr);
310   ierr = PetscObjectTypeCompare((PetscObject) B, MATSEQDENSE, &Bdense);CHKERRQ(ierr);
311   ierr = PetscObjectTypeCompare((PetscObject) C, MATSEQDENSE, &Cdense);CHKERRQ(ierr);
312   ierr = MatGetSize(B,&n_I,NULL);CHKERRQ(ierr);
313   if (n_I) {
314     if (!Bdense) {
315       ierr = MatConvert(B, MATSEQDENSE, MAT_INITIAL_MATRIX, &Bd);CHKERRQ(ierr);
316     } else {
317       Bd = B;
318     }
319 
320     if (isLU || isILU || isCHOL) {
321       Mat fact;
322       ierr = KSPSetUp(ksp);CHKERRQ(ierr);
323       ierr = PCFactorGetMatrix(pc, &fact);CHKERRQ(ierr);
324       ierr = MatDuplicate(Bd, MAT_DO_NOT_COPY_VALUES, &AinvBd);CHKERRQ(ierr);
325       ierr = MatMatSolve(fact, Bd, AinvBd);CHKERRQ(ierr);
326     } else {
327       PetscBool ex = PETSC_TRUE;
328 
329       if (ex) {
330         Mat Ainvd;
331 
332         ierr = PCComputeExplicitOperator(pc, &Ainvd);CHKERRQ(ierr);
333         ierr = MatMatMult(Ainvd, Bd, MAT_INITIAL_MATRIX, fill, &AinvBd);CHKERRQ(ierr);
334         ierr = MatDestroy(&Ainvd);CHKERRQ(ierr);
335       } else {
336         Vec         sol,rhs;
337         PetscScalar *arrayrhs,*arraysol;
338         PetscInt    i,nrhs,n;
339 
340         ierr = MatDuplicate(Bd, MAT_DO_NOT_COPY_VALUES, &AinvBd);CHKERRQ(ierr);
341         ierr = MatGetSize(Bd,&n,&nrhs);CHKERRQ(ierr);
342         ierr = MatDenseGetArray(Bd,&arrayrhs);CHKERRQ(ierr);
343         ierr = MatDenseGetArray(AinvBd,&arraysol);CHKERRQ(ierr);
344         ierr = KSPGetSolution(ksp,&sol);CHKERRQ(ierr);
345         ierr = KSPGetRhs(ksp,&rhs);CHKERRQ(ierr);
346         for (i=0;i<nrhs;i++) {
347           ierr = VecPlaceArray(rhs,arrayrhs+i*n);CHKERRQ(ierr);
348           ierr = VecPlaceArray(sol,arraysol+i*n);CHKERRQ(ierr);
349           ierr = KSPSolve(ksp,rhs,sol);CHKERRQ(ierr);
350           ierr = VecResetArray(rhs);CHKERRQ(ierr);
351           ierr = VecResetArray(sol);CHKERRQ(ierr);
352         }
353         ierr = MatDenseRestoreArray(Bd,&arrayrhs);CHKERRQ(ierr);
354         ierr = MatDenseRestoreArray(AinvBd,&arrayrhs);CHKERRQ(ierr);
355       }
356     }
357     if (!Bdense & !issym) {
358       ierr = MatDestroy(&Bd);CHKERRQ(ierr);
359     }
360 
361     if (!issym) {
362       if (!Cdense) {
363         ierr = MatConvert(C, MATSEQDENSE, MAT_INITIAL_MATRIX, &Cd);CHKERRQ(ierr);
364       } else {
365         Cd = C;
366       }
367       ierr = MatMatMult(Cd, AinvBd, reuse, fill, S);CHKERRQ(ierr);
368       if (!Cdense) {
369         ierr = MatDestroy(&Cd);CHKERRQ(ierr);
370       }
371     } else {
372       ierr = MatTransposeMatMult(Bd, AinvBd, reuse, fill, S);CHKERRQ(ierr);
373       if (!Bdense) {
374         ierr = MatDestroy(&Bd);CHKERRQ(ierr);
375       }
376     }
377     ierr = MatDestroy(&AinvBd);CHKERRQ(ierr);
378   }
379 
380   if (D) {
381     Mat       Dd;
382     PetscBool Ddense;
383 
384     ierr = PetscObjectTypeCompare((PetscObject)D,MATSEQDENSE,&Ddense);CHKERRQ(ierr);
385     if (!Ddense) {
386       ierr = MatConvert(D, MATSEQDENSE, MAT_INITIAL_MATRIX, &Dd);CHKERRQ(ierr);
387     } else {
388       Dd = D;
389     }
390     if (n_I) {
391       ierr = MatAYPX(*S,-1.0,Dd,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
392     } else {
393       if (reuse == MAT_INITIAL_MATRIX) {
394         ierr = MatDuplicate(Dd,MAT_COPY_VALUES,S);CHKERRQ(ierr);
395       } else {
396         ierr = MatCopy(Dd,*S,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
397       }
398     }
399     if (!Ddense) {
400       ierr = MatDestroy(&Dd);CHKERRQ(ierr);
401     }
402   } else {
403     ierr = MatScale(*S,-1.0);CHKERRQ(ierr);
404   }
405   PetscFunctionReturn(0);
406 }
407 
408 #undef __FUNCT__
409 #define __FUNCT__ "PCBDDCSubSchursSetUp"
410 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)
411 {
412   Mat                    F,A_II,A_IB,A_BI,A_BB,AE_II;
413   Mat                    S_all;
414   Mat                    global_schur_subsets,work_mat,*submats;
415   ISLocalToGlobalMapping l2gmap_subsets;
416   IS                     is_I,is_I_layer;
417   IS                     all_subsets,all_subsets_mult,all_subsets_n;
418   PetscInt               *nnz,*all_local_idx_N;
419   PetscInt               *auxnum1,*auxnum2;
420   PetscInt               i,subset_size,max_subset_size;
421   PetscInt               n_B,extra,local_size,global_size;
422   PetscBLASInt           B_N,B_ierr,B_lwork,*pivots;
423   PetscScalar            *Bwork;
424   PetscSubcomm           subcomm;
425   PetscMPIInt            color,rank;
426   MPI_Comm               comm_n;
427   PetscBool              deluxe = PETSC_TRUE, use_getr = PETSC_FALSE;
428   PetscErrorCode         ierr;
429 
430   PetscFunctionBegin;
431   /* update info in sub_schurs */
432   ierr = MatDestroy(&sub_schurs->A);CHKERRQ(ierr);
433   ierr = MatDestroy(&sub_schurs->S);CHKERRQ(ierr);
434   sub_schurs->is_hermitian = PETSC_FALSE;
435   sub_schurs->is_posdef = PETSC_FALSE;
436   if (Ain) {
437     PetscBool isseqaij;
438     /* determine if we are dealing with hermitian positive definite problems */
439 #if !defined(PETSC_USE_COMPLEX)
440     if (Ain->symmetric_set) {
441       sub_schurs->is_hermitian = Ain->symmetric;
442     }
443 #else
444     if (Ain->hermitian_set) {
445       sub_schurs->is_hermitian = Ain->hermitian;
446     }
447 #endif
448     if (Ain->spd_set) {
449       sub_schurs->is_posdef = Ain->spd;
450     }
451 
452     /* check */
453     ierr = PetscObjectTypeCompare((PetscObject)Ain,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
454     if (compute_Stilda && (!Ain->hermitian_set || !Ain->spd_set)) { /* these are lazy checks, maybe I should throw an error if they are not set */
455       PetscInt lsize;
456 
457       ierr = MatGetSize(Ain,&lsize,NULL);CHKERRQ(ierr);
458       if (lsize) {
459         PetscScalar val;
460         PetscReal   norm;
461         Vec         vec1,vec2,vec3;
462 
463         ierr = MatCreateVecs(Ain,&vec1,&vec2);CHKERRQ(ierr);
464         ierr = VecDuplicate(vec1,&vec3);CHKERRQ(ierr);
465         ierr = VecSetRandom(vec1,NULL);CHKERRQ(ierr);
466         ierr = MatMult(Ain,vec1,vec2);CHKERRQ(ierr);
467 #if !defined(PETSC_USE_COMPLEX)
468         ierr = MatMultTranspose(Ain,vec1,vec3);CHKERRQ(ierr);
469 #else
470         ierr = MatMultHermitianTranspose(Ain,vec1,vec3);CHKERRQ(ierr);
471 #endif
472         ierr = VecAXPY(vec3,-1.0,vec2);CHKERRQ(ierr);
473         ierr = VecNorm(vec3,NORM_INFINITY,&norm);CHKERRQ(ierr);
474         if (!Ain->hermitian_set) {
475           if (norm > PetscSqrtReal(PETSC_SMALL)) {
476             sub_schurs->is_hermitian = PETSC_FALSE;
477           } else {
478             sub_schurs->is_hermitian = PETSC_TRUE;
479           }
480         }
481         if (!Ain->spd_set && !benign_trick) {
482           ierr = VecDot(vec1,vec2,&val);CHKERRQ(ierr);
483           if (PetscRealPart(val) > 0. && PetscAbsReal(PetscImaginaryPart(val)) < PETSC_SMALL) sub_schurs->is_posdef = PETSC_TRUE;
484         }
485         ierr = VecDestroy(&vec1);CHKERRQ(ierr);
486         ierr = VecDestroy(&vec2);CHKERRQ(ierr);
487         ierr = VecDestroy(&vec3);CHKERRQ(ierr);
488       } else {
489         sub_schurs->is_hermitian = PETSC_TRUE;
490         sub_schurs->is_posdef = PETSC_TRUE;
491       }
492     }
493     if (isseqaij) {
494       ierr = PetscObjectReference((PetscObject)Ain);CHKERRQ(ierr);
495       sub_schurs->A = Ain;
496     } else {
497       ierr = MatConvert(Ain,MATSEQAIJ,MAT_INITIAL_MATRIX,&sub_schurs->A);CHKERRQ(ierr);
498     }
499   }
500 
501   ierr = PetscObjectReference((PetscObject)Sin);CHKERRQ(ierr);
502   sub_schurs->S = Sin;
503   if (sub_schurs->schur_explicit) {
504     sub_schurs->schur_explicit = (PetscBool)(!!sub_schurs->A);
505   }
506 
507   /* preliminary checks */
508   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");
509 
510   /* restrict work on active processes */
511   color = 0;
512   if (!sub_schurs->n_subs) color = 1; /* this can happen if we are in a multilevel case or if the subdomain is disconnected */
513   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&rank);CHKERRQ(ierr);
514   ierr = PetscSubcommCreate(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&subcomm);CHKERRQ(ierr);
515   ierr = PetscSubcommSetNumber(subcomm,2);CHKERRQ(ierr);
516   ierr = PetscSubcommSetTypeGeneral(subcomm,color,rank);CHKERRQ(ierr);
517   ierr = PetscCommDuplicate(PetscSubcommChild(subcomm),&comm_n,NULL);CHKERRQ(ierr);
518   ierr = PetscSubcommDestroy(&subcomm);CHKERRQ(ierr);
519   if (!sub_schurs->n_subs) {
520     ierr = PetscCommDestroy(&comm_n);CHKERRQ(ierr);
521     PetscFunctionReturn(0);
522   }
523   /* ierr = PetscCommDuplicate(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&comm_n,NULL);CHKERRQ(ierr); */
524 
525   /* get Schur complement matrices */
526   if (!sub_schurs->schur_explicit) {
527     Mat       tA_IB,tA_BI,tA_BB;
528     PetscBool isseqsbaij;
529     ierr = MatSchurComplementGetSubMatrices(sub_schurs->S,&A_II,NULL,&tA_IB,&tA_BI,&tA_BB);CHKERRQ(ierr);
530     ierr = PetscObjectTypeCompare((PetscObject)tA_BB,MATSEQSBAIJ,&isseqsbaij);CHKERRQ(ierr);
531     if (isseqsbaij) {
532       ierr = MatConvert(tA_BB,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_BB);CHKERRQ(ierr);
533       ierr = MatConvert(tA_IB,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_IB);CHKERRQ(ierr);
534       ierr = MatConvert(tA_BI,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_BI);CHKERRQ(ierr);
535     } else {
536       ierr = PetscObjectReference((PetscObject)tA_BB);CHKERRQ(ierr);
537       A_BB = tA_BB;
538       ierr = PetscObjectReference((PetscObject)tA_IB);CHKERRQ(ierr);
539       A_IB = tA_IB;
540       ierr = PetscObjectReference((PetscObject)tA_BI);CHKERRQ(ierr);
541       A_BI = tA_BI;
542     }
543   } else {
544     A_II = NULL;
545     A_IB = NULL;
546     A_BI = NULL;
547     A_BB = NULL;
548   }
549   S_all = NULL;
550 
551   /* determine interior problems */
552   ierr = ISGetLocalSize(sub_schurs->is_I,&i);CHKERRQ(ierr);
553   if (nlayers >= 0 && i) { /* Interior problems can be different from the original one */
554     PetscBT                touched;
555     const PetscInt*        idx_B;
556     PetscInt               n_I,n_B,n_local_dofs,n_prev_added,j,layer,*local_numbering;
557 
558     if (!xadj || !adjncy) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cannot request layering without adjacency");
559     /* get sizes */
560     ierr = ISGetLocalSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr);
561     ierr = ISGetLocalSize(sub_schurs->is_B,&n_B);CHKERRQ(ierr);
562 
563     ierr = PetscMalloc1(n_I+n_B,&local_numbering);CHKERRQ(ierr);
564     ierr = PetscBTCreate(n_I+n_B,&touched);CHKERRQ(ierr);
565     ierr = PetscBTMemzero(n_I+n_B,touched);CHKERRQ(ierr);
566 
567     /* all boundary dofs must be skipped when adding layers */
568     ierr = ISGetIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr);
569     for (j=0;j<n_B;j++) {
570       ierr = PetscBTSet(touched,idx_B[j]);CHKERRQ(ierr);
571     }
572     ierr = PetscMemcpy(local_numbering,idx_B,n_B*sizeof(PetscInt));CHKERRQ(ierr);
573     ierr = ISRestoreIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr);
574 
575     /* add prescribed number of layers of dofs */
576     n_local_dofs = n_B;
577     n_prev_added = n_B;
578     for (layer=0;layer<nlayers;layer++) {
579       PetscInt n_added;
580       if (n_local_dofs == n_I+n_B) break;
581       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);
582       ierr = PCBDDCAdjGetNextLayer_Private(local_numbering+n_local_dofs,n_prev_added,touched,xadj,adjncy,&n_added);CHKERRQ(ierr);
583       n_prev_added = n_added;
584       n_local_dofs += n_added;
585       if (!n_added) break;
586     }
587     ierr = PetscBTDestroy(&touched);CHKERRQ(ierr);
588 
589     /* IS for I layer dofs in original numbering */
590     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);
591     ierr = PetscFree(local_numbering);CHKERRQ(ierr);
592     ierr = ISSort(is_I_layer);CHKERRQ(ierr);
593     /* IS for I layer dofs in I numbering */
594     if (!sub_schurs->schur_explicit) {
595       ISLocalToGlobalMapping ItoNmap;
596       ierr = ISLocalToGlobalMappingCreateIS(sub_schurs->is_I,&ItoNmap);CHKERRQ(ierr);
597       ierr = ISGlobalToLocalMappingApplyIS(ItoNmap,IS_GTOLM_DROP,is_I_layer,&is_I);CHKERRQ(ierr);
598       ierr = ISLocalToGlobalMappingDestroy(&ItoNmap);CHKERRQ(ierr);
599 
600       /* II block */
601       ierr = MatGetSubMatrix(A_II,is_I,is_I,MAT_INITIAL_MATRIX,&AE_II);CHKERRQ(ierr);
602     }
603   } else {
604     PetscInt n_I;
605 
606     /* IS for I dofs in original numbering */
607     ierr = PetscObjectReference((PetscObject)sub_schurs->is_I);CHKERRQ(ierr);
608     is_I_layer = sub_schurs->is_I;
609 
610     /* IS for I dofs in I numbering (strided 1) */
611     if (!sub_schurs->schur_explicit) {
612       ierr = ISGetSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr);
613       ierr = ISCreateStride(PetscObjectComm((PetscObject)sub_schurs->is_I),n_I,0,1,&is_I);CHKERRQ(ierr);
614 
615       /* II block is the same */
616       ierr = PetscObjectReference((PetscObject)A_II);CHKERRQ(ierr);
617       AE_II = A_II;
618     }
619   }
620 
621   /* Get info on subset sizes and sum of all subsets sizes */
622   max_subset_size = 0;
623   local_size = 0;
624   for (i=0;i<sub_schurs->n_subs;i++) {
625     ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
626     max_subset_size = PetscMax(subset_size,max_subset_size);
627     local_size += subset_size;
628   }
629 
630   /* Work arrays for local indices */
631   extra = 0;
632   ierr = ISGetLocalSize(sub_schurs->is_B,&n_B);CHKERRQ(ierr);
633   if (sub_schurs->schur_explicit && is_I_layer) {
634     ierr = ISGetLocalSize(is_I_layer,&extra);CHKERRQ(ierr);
635   }
636   ierr = PetscMalloc1(n_B+extra,&all_local_idx_N);CHKERRQ(ierr);
637   if (extra) {
638     const PetscInt *idxs;
639     ierr = ISGetIndices(is_I_layer,&idxs);CHKERRQ(ierr);
640     ierr = PetscMemcpy(all_local_idx_N,idxs,extra*sizeof(PetscInt));CHKERRQ(ierr);
641     ierr = ISRestoreIndices(is_I_layer,&idxs);CHKERRQ(ierr);
642   }
643   ierr = PetscMalloc1(local_size,&nnz);CHKERRQ(ierr);
644   ierr = PetscMalloc1(sub_schurs->n_subs,&auxnum1);CHKERRQ(ierr);
645   ierr = PetscMalloc1(sub_schurs->n_subs,&auxnum2);CHKERRQ(ierr);
646 
647   /* Get local indices in local numbering */
648   local_size = 0;
649   for (i=0;i<sub_schurs->n_subs;i++) {
650     PetscInt j;
651     const    PetscInt *idxs;
652 
653     ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
654     ierr = ISGetIndices(sub_schurs->is_subs[i],&idxs);CHKERRQ(ierr);
655     /* start (smallest in global ordering) and multiplicity */
656     auxnum1[i] = idxs[0];
657     auxnum2[i] = subset_size;
658     /* subset indices in local numbering */
659     ierr = PetscMemcpy(all_local_idx_N+local_size+extra,idxs,subset_size*sizeof(PetscInt));CHKERRQ(ierr);
660     ierr = ISRestoreIndices(sub_schurs->is_subs[i],&idxs);CHKERRQ(ierr);
661     for (j=0;j<subset_size;j++) nnz[local_size+j] = subset_size;
662     local_size += subset_size;
663   }
664 
665   /* allocate extra workspace needed only for GETRI */
666   Bwork = NULL;
667   pivots = NULL;
668   if (local_size && !benign_trick && (!sub_schurs->is_hermitian || !sub_schurs->is_posdef)) {
669     PetscScalar lwork;
670 
671     use_getr = PETSC_TRUE;
672     B_lwork = -1;
673     ierr = PetscBLASIntCast(local_size,&B_N);CHKERRQ(ierr);
674     ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
675     PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,Bwork,&B_N,pivots,&lwork,&B_lwork,&B_ierr));
676     ierr = PetscFPTrapPop();CHKERRQ(ierr);
677     if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to GETRI Lapack routine %d",(int)B_ierr);
678     ierr = PetscBLASIntCast((PetscInt)PetscRealPart(lwork),&B_lwork);CHKERRQ(ierr);
679     ierr = PetscMalloc2(B_lwork,&Bwork,B_N,&pivots);CHKERRQ(ierr);
680   }
681 
682   /* prepare parallel matrices for summing up properly schurs on subsets */
683   ierr = ISCreateGeneral(comm_n,sub_schurs->n_subs,auxnum1,PETSC_OWN_POINTER,&all_subsets_n);CHKERRQ(ierr);
684   ierr = ISLocalToGlobalMappingApplyIS(sub_schurs->l2gmap,all_subsets_n,&all_subsets);CHKERRQ(ierr);
685   ierr = ISDestroy(&all_subsets_n);CHKERRQ(ierr);
686   ierr = ISCreateGeneral(comm_n,sub_schurs->n_subs,auxnum2,PETSC_OWN_POINTER,&all_subsets_mult);CHKERRQ(ierr);
687   ierr = PCBDDCSubsetNumbering(all_subsets,all_subsets_mult,&global_size,&all_subsets_n);CHKERRQ(ierr);
688   ierr = ISDestroy(&all_subsets);CHKERRQ(ierr);
689   ierr = ISDestroy(&all_subsets_mult);CHKERRQ(ierr);
690   ierr = ISGetLocalSize(all_subsets_n,&i);CHKERRQ(ierr);
691   if (i != local_size) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Invalid size of new subset! %D != %D",i,local_size);
692   ierr = ISLocalToGlobalMappingCreateIS(all_subsets_n,&l2gmap_subsets);CHKERRQ(ierr);
693   ierr = MatCreateIS(comm_n,1,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size,l2gmap_subsets,NULL,&work_mat);CHKERRQ(ierr);
694   ierr = ISLocalToGlobalMappingDestroy(&l2gmap_subsets);CHKERRQ(ierr);
695   ierr = MatCreate(PetscObjectComm((PetscObject)work_mat),&global_schur_subsets);CHKERRQ(ierr);
696   ierr = MatSetSizes(global_schur_subsets,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size);CHKERRQ(ierr);
697   ierr = MatSetType(global_schur_subsets,MATMPIAIJ);CHKERRQ(ierr);
698 
699   /* subset indices in local boundary numbering */
700   if (!sub_schurs->is_Ej_all) {
701     PetscInt *all_local_idx_B;
702 
703     ierr = PetscMalloc1(local_size,&all_local_idx_B);CHKERRQ(ierr);
704     ierr = ISGlobalToLocalMappingApply(sub_schurs->BtoNmap,IS_GTOLM_DROP,local_size,all_local_idx_N+extra,&subset_size,all_local_idx_B);CHKERRQ(ierr);
705     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);
706     ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size,all_local_idx_B,PETSC_OWN_POINTER,&sub_schurs->is_Ej_all);CHKERRQ(ierr);
707   }
708 
709   if (change) {
710     ISLocalToGlobalMapping BtoS;
711     IS                     change_primal_B;
712     IS                     change_primal_all;
713 
714     if (sub_schurs->change_primal_sub) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"This should not happen");
715     if (sub_schurs->change) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"This should not happen");
716     ierr = PetscMalloc1(sub_schurs->n_subs,&sub_schurs->change_primal_sub);CHKERRQ(ierr);
717     for (i=0;i<sub_schurs->n_subs;i++) {
718       ISLocalToGlobalMapping NtoS;
719       ierr = ISLocalToGlobalMappingCreateIS(sub_schurs->is_subs[i],&NtoS);CHKERRQ(ierr);
720       ierr = ISGlobalToLocalMappingApplyIS(NtoS,IS_GTOLM_DROP,change_primal,&sub_schurs->change_primal_sub[i]);CHKERRQ(ierr);
721       ierr = ISLocalToGlobalMappingDestroy(&NtoS);CHKERRQ(ierr);
722     }
723     ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,change_primal,&change_primal_B);CHKERRQ(ierr);
724     ierr = ISLocalToGlobalMappingCreateIS(sub_schurs->is_Ej_all,&BtoS);CHKERRQ(ierr);
725     ierr = ISGlobalToLocalMappingApplyIS(BtoS,IS_GTOLM_DROP,change_primal_B,&change_primal_all);CHKERRQ(ierr);
726     ierr = ISLocalToGlobalMappingDestroy(&BtoS);CHKERRQ(ierr);
727     ierr = ISDestroy(&change_primal_B);CHKERRQ(ierr);
728     ierr = PetscMalloc1(sub_schurs->n_subs,&sub_schurs->change);CHKERRQ(ierr);
729     for (i=0;i<sub_schurs->n_subs;i++) {
730       Mat change_sub;
731 
732       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
733       ierr = KSPCreate(PETSC_COMM_SELF,&sub_schurs->change[i]);CHKERRQ(ierr);
734       ierr = KSPSetType(sub_schurs->change[i],KSPPREONLY);CHKERRQ(ierr);
735       if (!sub_schurs->change_with_qr) {
736         ierr = MatGetSubMatrix(change,sub_schurs->is_subs[i],sub_schurs->is_subs[i],MAT_INITIAL_MATRIX,&change_sub);CHKERRQ(ierr);
737       } else {
738         Mat change_subt;
739         ierr = MatGetSubMatrix(change,sub_schurs->is_subs[i],sub_schurs->is_subs[i],MAT_INITIAL_MATRIX,&change_subt);CHKERRQ(ierr);
740         ierr = MatConvert(change_subt,MATSEQDENSE,MAT_INITIAL_MATRIX,&change_sub);CHKERRQ(ierr);
741         ierr = MatDestroy(&change_subt);CHKERRQ(ierr);
742       }
743       ierr = KSPSetOperators(sub_schurs->change[i],change_sub,change_sub);CHKERRQ(ierr);
744       ierr = MatDestroy(&change_sub);CHKERRQ(ierr);
745       ierr = KSPSetOptionsPrefix(sub_schurs->change[i],"sub_schurs_change_");CHKERRQ(ierr);
746     }
747     ierr = ISDestroy(&change_primal_all);CHKERRQ(ierr);
748   }
749 
750   /* Local matrix of all local Schur on subsets (transposed) */
751   if (!sub_schurs->S_Ej_all) {
752     ierr = MatCreate(PETSC_COMM_SELF,&sub_schurs->S_Ej_all);CHKERRQ(ierr);
753     ierr = MatSetSizes(sub_schurs->S_Ej_all,PETSC_DECIDE,PETSC_DECIDE,local_size,local_size);CHKERRQ(ierr);
754     ierr = MatSetType(sub_schurs->S_Ej_all,MATAIJ);CHKERRQ(ierr);
755     ierr = MatSeqAIJSetPreallocation(sub_schurs->S_Ej_all,0,nnz);CHKERRQ(ierr);
756   }
757 
758   /* Compute Schur complements explicitly */
759   F = NULL;
760   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 */
761     Mat         S_Ej_expl;
762     PetscScalar *work;
763     PetscInt    j,*dummy_idx;
764     PetscBool   Sdense;
765 
766     ierr = PetscMalloc2(max_subset_size,&dummy_idx,max_subset_size*max_subset_size,&work);CHKERRQ(ierr);
767     local_size = 0;
768     for (i=0;i<sub_schurs->n_subs;i++) {
769       IS  is_subset_B;
770       Mat AE_EE,AE_IE,AE_EI,S_Ej;
771 
772       /* subsets in original and boundary numbering */
773       ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,sub_schurs->is_subs[i],&is_subset_B);CHKERRQ(ierr);
774       /* EE block */
775       ierr = MatGetSubMatrix(A_BB,is_subset_B,is_subset_B,MAT_INITIAL_MATRIX,&AE_EE);CHKERRQ(ierr);
776       /* IE block */
777       ierr = MatGetSubMatrix(A_IB,is_I,is_subset_B,MAT_INITIAL_MATRIX,&AE_IE);CHKERRQ(ierr);
778       /* EI block */
779       if (sub_schurs->is_hermitian) {
780         ierr = MatCreateTranspose(AE_IE,&AE_EI);CHKERRQ(ierr);
781       } else {
782         ierr = MatGetSubMatrix(A_BI,is_subset_B,is_I,MAT_INITIAL_MATRIX,&AE_EI);CHKERRQ(ierr);
783       }
784       ierr = ISDestroy(&is_subset_B);CHKERRQ(ierr);
785       ierr = MatCreateSchurComplement(AE_II,AE_II,AE_IE,AE_EI,AE_EE,&S_Ej);CHKERRQ(ierr);
786       ierr = MatDestroy(&AE_EE);CHKERRQ(ierr);
787       ierr = MatDestroy(&AE_IE);CHKERRQ(ierr);
788       ierr = MatDestroy(&AE_EI);CHKERRQ(ierr);
789       if (AE_II == A_II) { /* we can reuse the same ksp */
790         KSP ksp;
791         ierr = MatSchurComplementGetKSP(sub_schurs->S,&ksp);CHKERRQ(ierr);
792         ierr = MatSchurComplementSetKSP(S_Ej,ksp);CHKERRQ(ierr);
793       } else { /* build new ksp object which inherits ksp and pc types from the original one */
794         KSP       origksp,schurksp;
795         PC        origpc,schurpc;
796         KSPType   ksp_type;
797         PetscInt  n_internal;
798         PetscBool ispcnone;
799 
800         ierr = MatSchurComplementGetKSP(sub_schurs->S,&origksp);CHKERRQ(ierr);
801         ierr = MatSchurComplementGetKSP(S_Ej,&schurksp);CHKERRQ(ierr);
802         ierr = KSPGetType(origksp,&ksp_type);CHKERRQ(ierr);
803         ierr = KSPSetType(schurksp,ksp_type);CHKERRQ(ierr);
804         ierr = KSPGetPC(schurksp,&schurpc);CHKERRQ(ierr);
805         ierr = KSPGetPC(origksp,&origpc);CHKERRQ(ierr);
806         ierr = PetscObjectTypeCompare((PetscObject)origpc,PCNONE,&ispcnone);CHKERRQ(ierr);
807         if (!ispcnone) {
808           PCType pc_type;
809           ierr = PCGetType(origpc,&pc_type);CHKERRQ(ierr);
810           ierr = PCSetType(schurpc,pc_type);CHKERRQ(ierr);
811         } else {
812           ierr = PCSetType(schurpc,PCLU);CHKERRQ(ierr);
813         }
814         ierr = ISGetSize(is_I,&n_internal);CHKERRQ(ierr);
815         if (n_internal) { /* UMFPACK gives error with 0 sized problems */
816           MatSolverPackage solver=NULL;
817           ierr = PCFactorGetMatSolverPackage(origpc,(const MatSolverPackage*)&solver);CHKERRQ(ierr);
818           if (solver) {
819             ierr = PCFactorSetMatSolverPackage(schurpc,solver);CHKERRQ(ierr);
820           }
821         }
822         ierr = KSPSetUp(schurksp);CHKERRQ(ierr);
823       }
824       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
825       ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&S_Ej_expl);CHKERRQ(ierr);
826       ierr = PCBDDCComputeExplicitSchur(S_Ej,sub_schurs->is_hermitian,MAT_REUSE_MATRIX,&S_Ej_expl);CHKERRQ(ierr);
827       ierr = PetscObjectTypeCompare((PetscObject)S_Ej_expl,MATSEQDENSE,&Sdense);CHKERRQ(ierr);
828       if (Sdense) {
829         for (j=0;j<subset_size;j++) {
830           dummy_idx[j]=local_size+j;
831         }
832         ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
833       } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not yet implemented for sparse matrices");
834       ierr = MatDestroy(&S_Ej);CHKERRQ(ierr);
835       ierr = MatDestroy(&S_Ej_expl);CHKERRQ(ierr);
836       local_size += subset_size;
837     }
838     ierr = PetscFree2(dummy_idx,work);CHKERRQ(ierr);
839     /* free */
840     ierr = ISDestroy(&is_I);CHKERRQ(ierr);
841     ierr = MatDestroy(&AE_II);CHKERRQ(ierr);
842     ierr = PetscFree(all_local_idx_N);CHKERRQ(ierr);
843   } else {
844     Mat         A,cs_AIB_mat = NULL,benign_AIIm1_ones_mat = NULL;
845     Vec         Dall = NULL;
846     IS          is_A_all,*is_p_r = NULL;
847     PetscScalar *work,*S_data,*schur_factor,infty = PETSC_MAX_REAL;
848     PetscInt    n,n_I,*dummy_idx,size_schur,size_active_schur,cum,cum2;
849     PetscBool   economic,solver_S,S_lower_triangular = PETSC_FALSE,factor_workaround;
850     char        solver[256];
851 
852     /* get sizes */
853     n_I = 0;
854     if (is_I_layer) {
855       ierr = ISGetLocalSize(is_I_layer,&n_I);CHKERRQ(ierr);
856     }
857     economic = PETSC_FALSE;
858     ierr = ISGetLocalSize(sub_schurs->is_I,&cum);CHKERRQ(ierr);
859     if (cum != n_I) economic = PETSC_TRUE;
860     ierr = MatGetLocalSize(sub_schurs->A,&n,NULL);CHKERRQ(ierr);
861     size_active_schur = local_size;
862 
863     /* import scaling vector */
864     if (scaling && compute_Stilda) {
865       const PetscScalar *array;
866       PetscScalar       *array2;
867       const PetscInt    *idxs;
868       PetscInt          i;
869 
870       ierr = ISGetIndices(sub_schurs->is_Ej_all,&idxs);CHKERRQ(ierr);
871       ierr = VecCreateSeq(PETSC_COMM_SELF,size_active_schur,&Dall);CHKERRQ(ierr);
872       ierr = VecGetArrayRead(scaling,&array);CHKERRQ(ierr);
873       ierr = VecGetArray(Dall,&array2);CHKERRQ(ierr);
874       for (i=0;i<size_active_schur;i++) array2[i] = array[idxs[i]];
875       ierr = VecRestoreArray(Dall,&array2);CHKERRQ(ierr);
876       ierr = VecRestoreArrayRead(scaling,&array);CHKERRQ(ierr);
877       ierr = ISRestoreIndices(sub_schurs->is_Ej_all,&idxs);CHKERRQ(ierr);
878       deluxe = PETSC_FALSE;
879     }
880 
881     /* size active schurs does not count any dirichlet or vertex dof on the interface */
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(NULL,"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,Afake;
1096       Vec                vec1_B;
1097       PCBDDCReuseSolvers msolv_ctx;
1098       PetscInt           n_R;
1099 
1100       if (sub_schurs->reuse_solver) {
1101         ierr = PCBDDCReuseSolversReset(sub_schurs->reuse_solver);CHKERRQ(ierr);
1102       } else {
1103         ierr = PetscNew(&sub_schurs->reuse_solver);CHKERRQ(ierr);
1104       }
1105       msolv_ctx = sub_schurs->reuse_solver;
1106       ierr = MatSchurComplementGetSubMatrices(sub_schurs->S,&A_II,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
1107       ierr = PetscObjectReference((PetscObject)F);CHKERRQ(ierr);
1108       msolv_ctx->F = F;
1109       ierr = MatCreateVecs(F,&msolv_ctx->sol,NULL);CHKERRQ(ierr);
1110       /* currently PETSc has no support for MatSolve(F,x,x), so cheat and let rhs and sol share the same memory */
1111       {
1112         PetscScalar *array;
1113         PetscInt    n;
1114 
1115         ierr = VecGetLocalSize(msolv_ctx->sol,&n);CHKERRQ(ierr);
1116         ierr = VecGetArray(msolv_ctx->sol,&array);CHKERRQ(ierr);
1117         ierr = VecCreateSeqWithArray(PetscObjectComm((PetscObject)msolv_ctx->sol),1,n,array,&msolv_ctx->rhs);CHKERRQ(ierr);
1118         ierr = VecRestoreArray(msolv_ctx->sol,&array);CHKERRQ(ierr);
1119       }
1120       msolv_ctx->has_vertices = factor_workaround;
1121 
1122       /* interior solver */
1123       ierr = PCCreate(PetscObjectComm((PetscObject)A_II),&msolv_ctx->interior_solver);CHKERRQ(ierr);
1124       ierr = PCSetOperators(msolv_ctx->interior_solver,A_II,A_II);CHKERRQ(ierr);
1125       ierr = PCSetType(msolv_ctx->interior_solver,PCSHELL);CHKERRQ(ierr);
1126       ierr = PCShellSetContext(msolv_ctx->interior_solver,msolv_ctx);CHKERRQ(ierr);
1127       ierr = PCShellSetApply(msolv_ctx->interior_solver,PCBDDCReuseSolvers_Interior);CHKERRQ(ierr);
1128       ierr = PCShellSetApplyTranspose(msolv_ctx->interior_solver,PCBDDCReuseSolvers_InteriorTranspose);CHKERRQ(ierr);
1129 
1130       /* correction solver */
1131       ierr = PCCreate(PetscObjectComm((PetscObject)A_II),&msolv_ctx->correction_solver);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 = ISGetLocalSize(msolv_ctx->is_R,&n_R);CHKERRQ(ierr);
1164       ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,n_R,n_R,0,NULL,&Afake);CHKERRQ(ierr);
1165       ierr = PCSetOperators(msolv_ctx->correction_solver,Afake,Afake);CHKERRQ(ierr);
1166       ierr = MatDestroy(&Afake);CHKERRQ(ierr);
1167       ierr = VecDestroy(&vec1_B);CHKERRQ(ierr);
1168 
1169       /* communicate benign info to solver context */
1170       if (benign_n) {
1171         PetscScalar *array;
1172 
1173         msolv_ctx->benign_n = benign_n;
1174         msolv_ctx->benign_zerodiag_subs = is_p_r;
1175         ierr = PetscMalloc1(benign_n,&msolv_ctx->benign_save_vals);CHKERRQ(ierr);
1176         msolv_ctx->benign_csAIB = cs_AIB_mat;
1177         ierr = MatCreateVecs(cs_AIB_mat,NULL,&msolv_ctx->benign_corr_work);CHKERRQ(ierr);
1178         ierr = VecGetArray(msolv_ctx->benign_corr_work,&array);CHKERRQ(ierr);
1179         ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,size_schur,array,&msolv_ctx->benign_dummy_schur_vec);CHKERRQ(ierr);
1180         ierr = VecRestoreArray(msolv_ctx->benign_corr_work,&array);CHKERRQ(ierr);
1181         msolv_ctx->benign_AIIm1ones = benign_AIIm1_ones_mat;
1182       }
1183     }
1184     ierr = MatDestroy(&A);CHKERRQ(ierr);
1185     ierr = ISDestroy(&is_A_all);CHKERRQ(ierr);
1186 
1187     /* Work arrays */
1188     ierr = PetscMalloc2(max_subset_size,&dummy_idx,max_subset_size*max_subset_size,&work);CHKERRQ(ierr);
1189 
1190     /* matrices for adaptive selection */
1191     if (compute_Stilda) {
1192       if (!sub_schurs->sum_S_Ej_tilda_all) {
1193         ierr = MatCreate(PETSC_COMM_SELF,&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
1194         ierr = MatSetSizes(sub_schurs->sum_S_Ej_tilda_all,PETSC_DECIDE,PETSC_DECIDE,size_active_schur,size_active_schur);CHKERRQ(ierr);
1195         ierr = MatSetType(sub_schurs->sum_S_Ej_tilda_all,MATAIJ);CHKERRQ(ierr);
1196         ierr = MatSeqAIJSetPreallocation(sub_schurs->sum_S_Ej_tilda_all,0,nnz);CHKERRQ(ierr);
1197       }
1198       if (!sub_schurs->sum_S_Ej_inv_all && deluxe) {
1199         ierr = MatCreate(PETSC_COMM_SELF,&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
1200         ierr = MatSetSizes(sub_schurs->sum_S_Ej_inv_all,PETSC_DECIDE,PETSC_DECIDE,size_active_schur,size_active_schur);CHKERRQ(ierr);
1201         ierr = MatSetType(sub_schurs->sum_S_Ej_inv_all,MATAIJ);CHKERRQ(ierr);
1202         ierr = MatSeqAIJSetPreallocation(sub_schurs->sum_S_Ej_inv_all,0,nnz);CHKERRQ(ierr);
1203       }
1204     }
1205 
1206     /* S_Ej_all */
1207     cum = cum2 = 0;
1208     ierr = MatDenseGetArray(S_all,&S_data);CHKERRQ(ierr);
1209     for (i=0;i<sub_schurs->n_subs;i++) {
1210       PetscInt j;
1211 
1212       /* get S_E */
1213       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
1214       if (S_lower_triangular) { /* I need to expand the upper triangular data (column oriented) */
1215         PetscInt k;
1216         for (k=0;k<subset_size;k++) {
1217           for (j=k;j<subset_size;j++) {
1218             work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
1219             work[j*subset_size+k] = PetscConj(S_data[cum2+k*size_schur+j]);
1220           }
1221         }
1222       } else { /* just copy to workspace */
1223         PetscInt k;
1224         for (k=0;k<subset_size;k++) {
1225           for (j=0;j<subset_size;j++) {
1226             work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
1227           }
1228         }
1229       }
1230       /* insert S_E values */
1231       for (j=0;j<subset_size;j++) dummy_idx[j] = cum+j;
1232       if (sub_schurs->change) {
1233         Mat            change_sub,SEj,T;
1234         const PetscInt *idxs;
1235         PetscInt       nz;
1236 
1237         /* change basis */
1238         ierr = KSPGetOperators(sub_schurs->change[i],&change_sub,NULL);CHKERRQ(ierr);
1239         ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&SEj);CHKERRQ(ierr);
1240         if (!sub_schurs->change_with_qr) { /* currently there's no support for PtAP with P SeqAIJ */
1241           Mat T2;
1242           ierr = MatTransposeMatMult(change_sub,SEj,MAT_INITIAL_MATRIX,1.0,&T2);CHKERRQ(ierr);
1243           ierr = MatMatMult(T2,change_sub,MAT_INITIAL_MATRIX,1.0,&T);CHKERRQ(ierr);
1244           ierr = MatConvert(T,MATSEQDENSE,MAT_INPLACE_MATRIX,&T);CHKERRQ(ierr);
1245           ierr = MatDestroy(&T2);CHKERRQ(ierr);
1246         } else {
1247           ierr = MatPtAP(SEj,change_sub,MAT_INITIAL_MATRIX,1.0,&T);CHKERRQ(ierr);
1248         }
1249         ierr = MatCopy(T,SEj,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
1250         ierr = MatDestroy(&T);CHKERRQ(ierr);
1251         /* currently there's no support for ZeroRowsColumns with A dense */
1252         ierr = ISGetIndices(sub_schurs->change_primal_sub[i],&idxs);CHKERRQ(ierr);
1253         ierr = ISGetLocalSize(sub_schurs->change_primal_sub[i],&nz);CHKERRQ(ierr);
1254         for (j=0;j<nz;j++) {
1255           PetscInt k;
1256 
1257           ierr = PetscMemzero(work+subset_size*idxs[j],subset_size*sizeof(PetscScalar));CHKERRQ(ierr);
1258           for (k=0;k<subset_size;k++) work[idxs[j]+k*subset_size] = 0.;
1259           work[idxs[j]*(subset_size+1)] = 1.0;
1260         }
1261         ierr = ISRestoreIndices(sub_schurs->change_primal_sub[i],&idxs);CHKERRQ(ierr);
1262         ierr = MatDestroy(&SEj);CHKERRQ(ierr);
1263       }
1264       if (deluxe) {
1265         ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
1266         /* if adaptivity is requested, invert S_E blocks */
1267         if (compute_Stilda) {
1268           ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr);
1269           ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
1270           if (!use_getr) { /* TODO add sytrf/i for symmetric non hermitian */
1271             PetscInt k;
1272             PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,work,&B_N,&B_ierr));
1273             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
1274             PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,work,&B_N,&B_ierr));
1275             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
1276             for (k=0;k<subset_size;k++) {
1277               for (j=k;j<subset_size;j++) {
1278                 work[j*subset_size+k] = work[k*subset_size+j];
1279               }
1280             }
1281           } else {
1282             PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,work,&B_N,pivots,&B_ierr));
1283             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
1284             PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,work,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
1285             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
1286           }
1287           ierr = PetscFPTrapPop();CHKERRQ(ierr);
1288           ierr = MatSetValues(sub_schurs->sum_S_Ej_inv_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
1289         }
1290       } else if (compute_Stilda) { /* not using deluxe */
1291         Mat         SEj;
1292         Vec         D;
1293         PetscScalar *array;
1294 
1295         ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&SEj);CHKERRQ(ierr);
1296         ierr = VecGetArray(Dall,&array);CHKERRQ(ierr);
1297         ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,subset_size,array+cum,&D);CHKERRQ(ierr);
1298         ierr = VecRestoreArray(Dall,&array);CHKERRQ(ierr);
1299         ierr = MatDiagonalScale(SEj,D,D);CHKERRQ(ierr);
1300         ierr = MatDestroy(&SEj);CHKERRQ(ierr);
1301         ierr = VecDestroy(&D);CHKERRQ(ierr);
1302         ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
1303       }
1304       cum += subset_size;
1305       cum2 += subset_size*(size_schur + 1);
1306     }
1307     ierr = MatDenseRestoreArray(S_all,&S_data);CHKERRQ(ierr);
1308 
1309     if (solver_S) {
1310       ierr = MatFactorRestoreSchurComplement(F,&S_all);CHKERRQ(ierr);
1311     }
1312 
1313     schur_factor = NULL;
1314     if (compute_Stilda && size_active_schur) {
1315 
1316       if (sub_schurs->n_subs == 1 && size_schur == size_active_schur && deluxe) { /* we already computed the inverse */
1317         PetscInt j;
1318         for (j=0;j<size_schur;j++) dummy_idx[j] = j;
1319         ierr = MatSetValues(sub_schurs->sum_S_Ej_tilda_all,size_schur,dummy_idx,size_schur,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
1320       } else {
1321         Mat S_all_inv=NULL;
1322         if (solver_S) { /* use MatFactor calls to invert S */
1323             /* let MatFactor factor the Schur complement */
1324           ierr = MatFactorFactorizeSchurComplement(F);CHKERRQ(ierr);
1325           /* for adaptive selection we need S^-1; for solver reusage we need S_\Delta\Delta^-1.
1326              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 */
1327           if (factor_workaround) {
1328             PetscScalar *data;
1329             PetscInt nd = 0;
1330 
1331             if (!sub_schurs->is_posdef) {
1332               SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor update not yet implemented for non SPD matrices");
1333             }
1334             ierr = MatFactorGetSchurComplement(F,&S_all_inv);CHKERRQ(ierr);
1335             ierr = MatDenseGetArray(S_all_inv,&data);CHKERRQ(ierr);
1336             if (sub_schurs->is_dir) { /* dirichlet dofs could have different scalings */
1337               ierr = ISGetLocalSize(sub_schurs->is_dir,&nd);CHKERRQ(ierr);
1338             }
1339             ierr = PetscMalloc1((size_active_schur*(size_active_schur+1))/2+nd,&schur_factor);CHKERRQ(ierr);
1340             cum = 0;
1341             for (i=0;i<size_active_schur;i++) {
1342               ierr = PetscMemcpy(schur_factor+cum,data+i*(size_schur+1),(size_active_schur-i)*sizeof(PetscScalar));CHKERRQ(ierr);
1343               cum += size_active_schur-i;
1344             }
1345             for (i=0;i<nd;i++) schur_factor[cum+i] = data[(i+size_active_schur)*(size_schur+1)];
1346             /* invert without calling MatFactorInvertSchurComplement, since we are hacking */
1347             ierr = PetscBLASIntCast(size_schur,&B_N);CHKERRQ(ierr);
1348             ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
1349             PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,data,&B_N,&B_ierr));
1350             ierr = PetscFPTrapPop();CHKERRQ(ierr);
1351             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
1352             ierr = MatDenseRestoreArray(S_all_inv,&data);CHKERRQ(ierr);
1353           } else {
1354             ierr = MatFactorInvertSchurComplement(F);CHKERRQ(ierr);
1355             ierr = MatFactorGetSchurComplement(F,&S_all_inv);CHKERRQ(ierr);
1356           }
1357         } else { /* we need to invert explicitly since we are not using MatFactor for S */
1358           ierr = PetscObjectReference((PetscObject)S_all);CHKERRQ(ierr);
1359           S_all_inv = S_all;
1360           ierr = MatDenseGetArray(S_all_inv,&S_data);CHKERRQ(ierr);
1361           ierr = PetscBLASIntCast(size_schur,&B_N);CHKERRQ(ierr);
1362           ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
1363           if (!use_getr) { /* TODO add sytrf/i for symmetric non hermitian */
1364             PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,S_data,&B_N,&B_ierr));
1365             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
1366             PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,S_data,&B_N,&B_ierr));
1367             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
1368           } else {
1369             PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,S_data,&B_N,pivots,&B_ierr));
1370             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
1371             PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,S_data,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
1372             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
1373           }
1374           ierr = PetscFPTrapPop();CHKERRQ(ierr);
1375           ierr = MatDenseRestoreArray(S_all_inv,&S_data);CHKERRQ(ierr);
1376         }
1377         /* S_Ej_tilda_all */
1378         cum = cum2 = 0;
1379         ierr = MatDenseGetArray(S_all_inv,&S_data);CHKERRQ(ierr);
1380         for (i=0;i<sub_schurs->n_subs;i++) {
1381           PetscInt j;
1382 
1383           ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
1384           /* get (St^-1)_E */
1385           /* Unless we are changing the variables, I don't need to expand to upper triangular since St^-1
1386              will be properly accessed later during adaptive selection */
1387           if (S_lower_triangular) {
1388             PetscInt k;
1389             if (sub_schurs->change) {
1390               for (k=0;k<subset_size;k++) {
1391                 for (j=k;j<subset_size;j++) {
1392                   work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
1393                   work[j*subset_size+k] = work[k*subset_size+j];
1394                 }
1395               }
1396             } else {
1397               for (k=0;k<subset_size;k++) {
1398                 for (j=k;j<subset_size;j++) {
1399                   work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
1400                 }
1401               }
1402             }
1403           } else {
1404             PetscInt k;
1405             for (k=0;k<subset_size;k++) {
1406               for (j=0;j<subset_size;j++) {
1407                 work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
1408               }
1409             }
1410           }
1411           if (sub_schurs->change) {
1412             Mat change_sub,SEj,T;
1413 
1414             /* change basis */
1415             ierr = KSPGetOperators(sub_schurs->change[i],&change_sub,NULL);CHKERRQ(ierr);
1416             ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&SEj);CHKERRQ(ierr);
1417             if (!sub_schurs->change_with_qr) { /* currently there's no support for PtAP with P SeqAIJ */
1418               Mat T2;
1419               ierr = MatTransposeMatMult(change_sub,SEj,MAT_INITIAL_MATRIX,1.0,&T2);CHKERRQ(ierr);
1420               ierr = MatMatMult(T2,change_sub,MAT_INITIAL_MATRIX,1.0,&T);CHKERRQ(ierr);
1421               ierr = MatDestroy(&T2);CHKERRQ(ierr);
1422               ierr = MatConvert(T,MATSEQDENSE,MAT_INPLACE_MATRIX,&T);CHKERRQ(ierr);
1423             } else {
1424               ierr = MatPtAP(SEj,change_sub,MAT_INITIAL_MATRIX,1.0,&T);CHKERRQ(ierr);
1425             }
1426             ierr = MatCopy(T,SEj,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
1427             ierr = MatDestroy(&T);CHKERRQ(ierr);
1428             /* set diagonal entry to a very large value to pick the basis we are eliminating as the first eigenvectors with adaptive selection */
1429             {
1430               const PetscInt *idxs;
1431               PetscInt n,j;
1432               ierr = ISGetLocalSize(sub_schurs->change_primal_sub[i],&n);CHKERRQ(ierr);
1433               ierr = ISGetIndices(sub_schurs->change_primal_sub[i],&idxs);CHKERRQ(ierr);
1434               for (j=0;j<n;j++) {
1435                 PetscInt k;
1436                 for (k=0;k<subset_size;k++) {
1437                   work[k+idxs[j]*subset_size] = work[idxs[j]+k*subset_size] = 0.;
1438                 }
1439                 work[idxs[j] + idxs[j]*subset_size] = 1./PETSC_SMALL;
1440               }
1441               ierr = ISRestoreIndices(sub_schurs->change_primal_sub[i],&idxs);CHKERRQ(ierr);
1442             }
1443             ierr = MatDestroy(&SEj);CHKERRQ(ierr);
1444           }
1445           for (j=0;j<subset_size;j++) dummy_idx[j] = cum+j;
1446           ierr = MatSetValues(sub_schurs->sum_S_Ej_tilda_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
1447           cum += subset_size;
1448           cum2 += subset_size*(size_schur + 1);
1449         }
1450         ierr = MatDenseRestoreArray(S_all_inv,&S_data);CHKERRQ(ierr);
1451         if (solver_S) {
1452           ierr = MatFactorRestoreSchurComplement(F,&S_all_inv);CHKERRQ(ierr);
1453         }
1454         ierr = MatDestroy(&S_all_inv);CHKERRQ(ierr);
1455       }
1456 
1457       /* move back factors to Schur data into F */
1458       if (factor_workaround) {
1459         Mat S_tmp;
1460         PetscInt nv,nd=0;
1461 
1462         if (!solver_S) {
1463           SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"This should not happen");
1464         }
1465         ierr = ISGetLocalSize(sub_schurs->is_vertices,&nv);CHKERRQ(ierr);
1466         ierr = MatFactorGetSchurComplement(F,&S_tmp);CHKERRQ(ierr);
1467         if (sub_schurs->is_posdef) {
1468           PetscScalar *data;
1469 
1470           ierr = MatZeroEntries(S_tmp);CHKERRQ(ierr);
1471           ierr = MatDenseGetArray(S_tmp,&data);CHKERRQ(ierr);
1472 
1473           if (S_lower_triangular) {
1474             cum = 0;
1475             for (i=0;i<size_active_schur;i++) {
1476               ierr = PetscMemcpy(data+i*(size_schur+1),schur_factor+cum,(size_active_schur-i)*sizeof(PetscScalar));CHKERRQ(ierr);
1477               cum += size_active_schur-i;
1478 	    }
1479           } else {
1480             ierr = PetscMemcpy(data,schur_factor,size_schur*size_schur*sizeof(PetscScalar));CHKERRQ(ierr);
1481           }
1482           if (sub_schurs->is_dir) {
1483             ierr = ISGetLocalSize(sub_schurs->is_dir,&nd);CHKERRQ(ierr);
1484 	    for (i=0;i<nd;i++) {
1485 	      data[(i+size_active_schur)*(size_schur+1)] = schur_factor[cum+i];
1486 	    }
1487           }
1488           /* workaround: since I cannot modify the matrices used inside the solvers for the forward and backward substitutions,
1489              set the diagonal entry of the Schur factor to a very large value */
1490           for (i=size_active_schur+nd;i<size_schur;i++) {
1491             data[i*(size_schur+1)] = infty;
1492           }
1493           ierr = MatDenseRestoreArray(S_tmp,&data);CHKERRQ(ierr);
1494         } else {
1495           SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor update not yet implemented for non SPD matrices");
1496         }
1497         ierr = MatFactorRestoreSchurComplement(F,&S_tmp);CHKERRQ(ierr);
1498       }
1499     } else if (factor_workaround) { /* we need to eliminate any unneeded coupling */
1500       PetscScalar *data;
1501       PetscInt    nd = 0;
1502 
1503       if (sub_schurs->is_dir) { /* dirichlet dofs could have different scalings */
1504         ierr = ISGetLocalSize(sub_schurs->is_dir,&nd);CHKERRQ(ierr);
1505       }
1506       ierr = MatFactorGetSchurComplement(F,&S_all);CHKERRQ(ierr);
1507       ierr = MatDenseGetArray(S_all,&data);CHKERRQ(ierr);
1508       for (i=0;i<size_active_schur;i++) {
1509         ierr = PetscMemzero(data+i*size_schur+size_active_schur,(size_schur-size_active_schur)*sizeof(PetscScalar));CHKERRQ(ierr);
1510       }
1511       for (i=size_active_schur+nd;i<size_schur;i++) {
1512         ierr = PetscMemzero(data+i*size_schur+size_active_schur,(size_schur-size_active_schur)*sizeof(PetscScalar));CHKERRQ(ierr);
1513         data[i*(size_schur+1)] = infty;
1514       }
1515       ierr = MatDenseRestoreArray(S_all,&data);CHKERRQ(ierr);
1516       ierr = MatFactorRestoreSchurComplement(F,&S_all);CHKERRQ(ierr);
1517     }
1518     ierr = PetscFree2(dummy_idx,work);CHKERRQ(ierr);
1519     ierr = PetscFree(schur_factor);CHKERRQ(ierr);
1520     ierr = VecDestroy(&Dall);CHKERRQ(ierr);
1521   }
1522   ierr = PetscFree(nnz);CHKERRQ(ierr);
1523   ierr = MatDestroy(&F);CHKERRQ(ierr);
1524   ierr = ISDestroy(&is_I_layer);CHKERRQ(ierr);
1525   ierr = MatDestroy(&S_all);CHKERRQ(ierr);
1526   ierr = MatDestroy(&A_BB);CHKERRQ(ierr);
1527   ierr = MatDestroy(&A_IB);CHKERRQ(ierr);
1528   ierr = MatDestroy(&A_BI);CHKERRQ(ierr);
1529   ierr = MatAssemblyBegin(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1530   ierr = MatAssemblyEnd(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1531   if (compute_Stilda) {
1532     ierr = MatAssemblyBegin(sub_schurs->sum_S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1533     ierr = MatAssemblyEnd(sub_schurs->sum_S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1534     if (deluxe) {
1535       ierr = MatAssemblyBegin(sub_schurs->sum_S_Ej_inv_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1536       ierr = MatAssemblyEnd(sub_schurs->sum_S_Ej_inv_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1537     }
1538   }
1539 
1540   /* Global matrix of all assembled Schur on subsets */
1541   ierr = MatISSetLocalMat(work_mat,sub_schurs->S_Ej_all);CHKERRQ(ierr);
1542   ierr = MatISSetMPIXAIJPreallocation_Private(work_mat,global_schur_subsets,PETSC_TRUE);CHKERRQ(ierr);
1543   ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
1544 
1545   /* Get local part of (\sum_j S_Ej) */
1546   if (!sub_schurs->sum_S_Ej_all) {
1547     ierr = MatGetSubMatrices(global_schur_subsets,1,&all_subsets_n,&all_subsets_n,MAT_INITIAL_MATRIX,&submats);CHKERRQ(ierr);
1548     sub_schurs->sum_S_Ej_all = submats[0];
1549   } else {
1550     ierr = PetscMalloc1(1,&submats);CHKERRQ(ierr);
1551     submats[0] = sub_schurs->sum_S_Ej_all;
1552     ierr = MatGetSubMatrices(global_schur_subsets,1,&all_subsets_n,&all_subsets_n,MAT_REUSE_MATRIX,&submats);CHKERRQ(ierr);
1553   }
1554 
1555   /* Get local part of (\sum_j S^-1_Ej) (\sum_j St^-1_Ej) */
1556   if (compute_Stilda) {
1557     ierr = MatISSetLocalMat(work_mat,sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
1558     ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
1559     submats[0] = sub_schurs->sum_S_Ej_tilda_all;
1560     ierr = MatGetSubMatrices(global_schur_subsets,1,&all_subsets_n,&all_subsets_n,MAT_REUSE_MATRIX,&submats);CHKERRQ(ierr);
1561     if (deluxe) {
1562       ierr = MatISSetLocalMat(work_mat,sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
1563       ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
1564       submats[0] = sub_schurs->sum_S_Ej_inv_all;
1565       ierr = MatGetSubMatrices(global_schur_subsets,1,&all_subsets_n,&all_subsets_n,MAT_REUSE_MATRIX,&submats);CHKERRQ(ierr);
1566     } else {
1567       PetscScalar *array;
1568       PetscInt    cum;
1569 
1570       ierr = MatSeqAIJGetArray(sub_schurs->sum_S_Ej_tilda_all,&array);CHKERRQ(ierr);
1571       cum = 0;
1572       for (i=0;i<sub_schurs->n_subs;i++) {
1573         ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
1574         ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr);
1575         ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
1576         PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,array+cum,&B_N,&B_ierr));
1577         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
1578         PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,array+cum,&B_N,&B_ierr));
1579         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
1580         ierr = PetscFPTrapPop();CHKERRQ(ierr);
1581         cum += subset_size*subset_size;
1582       }
1583       ierr = MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_tilda_all,&array);CHKERRQ(ierr);
1584       ierr = PetscObjectReference((PetscObject)sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
1585       ierr = MatDestroy(&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
1586       sub_schurs->sum_S_Ej_inv_all = sub_schurs->sum_S_Ej_all;
1587     }
1588   }
1589 
1590   /* free workspace */
1591   ierr = PetscFree(submats);CHKERRQ(ierr);
1592   ierr = PetscFree2(Bwork,pivots);CHKERRQ(ierr);
1593   ierr = MatDestroy(&global_schur_subsets);CHKERRQ(ierr);
1594   ierr = MatDestroy(&work_mat);CHKERRQ(ierr);
1595   ierr = ISDestroy(&all_subsets_n);CHKERRQ(ierr);
1596   ierr = PetscCommDestroy(&comm_n);CHKERRQ(ierr);
1597   PetscFunctionReturn(0);
1598 }
1599 
1600 #undef __FUNCT__
1601 #define __FUNCT__ "PCBDDCSubSchursInit"
1602 PetscErrorCode PCBDDCSubSchursInit(PCBDDCSubSchurs sub_schurs, IS is_I, IS is_B, PCBDDCGraph graph, ISLocalToGlobalMapping BtoNmap)
1603 {
1604   IS              *faces,*edges,*all_cc,vertices;
1605   PetscInt        i,n_faces,n_edges,n_all_cc;
1606   PetscBool       is_sorted;
1607   PetscErrorCode  ierr;
1608 
1609   PetscFunctionBegin;
1610   ierr = ISSorted(is_I,&is_sorted);CHKERRQ(ierr);
1611   if (!is_sorted) SETERRQ(PetscObjectComm((PetscObject)is_I),PETSC_ERR_PLIB,"IS for I dofs should be shorted");
1612   ierr = ISSorted(is_B,&is_sorted);CHKERRQ(ierr);
1613   if (!is_sorted) SETERRQ(PetscObjectComm((PetscObject)is_B),PETSC_ERR_PLIB,"IS for B dofs should be shorted");
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       ierr = ISDestroy(&sub_schurs->change_primal_sub[i]);CHKERRQ(ierr);
1718     }
1719   }
1720   ierr = PetscFree(sub_schurs->change);CHKERRQ(ierr);
1721   ierr = PetscFree(sub_schurs->change_primal_sub);CHKERRQ(ierr);
1722   sub_schurs->n_subs = 0;
1723   PetscFunctionReturn(0);
1724 }
1725 
1726 #undef __FUNCT__
1727 #define __FUNCT__ "PCBDDCAdjGetNextLayer_Private"
1728 PETSC_STATIC_INLINE PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt* queue_tip,PetscInt n_prev,PetscBT touched,PetscInt* xadj,PetscInt* adjncy,PetscInt* n_added)
1729 {
1730   PetscInt       i,j,n;
1731   PetscErrorCode ierr;
1732 
1733   PetscFunctionBegin;
1734   n = 0;
1735   for (i=-n_prev;i<0;i++) {
1736     PetscInt start_dof = queue_tip[i];
1737     for (j=xadj[start_dof];j<xadj[start_dof+1];j++) {
1738       PetscInt dof = adjncy[j];
1739       if (!PetscBTLookup(touched,dof)) {
1740         ierr = PetscBTSet(touched,dof);CHKERRQ(ierr);
1741         queue_tip[n] = dof;
1742         n++;
1743       }
1744     }
1745   }
1746   *n_added = n;
1747   PetscFunctionReturn(0);
1748 }
1749