xref: /petsc/src/ksp/pc/impls/bddc/bddcschurs.c (revision 9d54b7f467170aa44380999bc2040dbb24de4386)
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, PetscBool faster_deluxe, 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,*change_primal_sub;
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   change_primal_sub = NULL;
710   if (change) {
711     ISLocalToGlobalMapping BtoS;
712     IS                     change_primal_B;
713     IS                     change_primal_all;
714 
715     ierr = PetscMalloc1(sub_schurs->n_subs,&change_primal_sub);CHKERRQ(ierr);
716     for (i=0;i<sub_schurs->n_subs;i++) {
717       ISLocalToGlobalMapping NtoS;
718       ierr = ISLocalToGlobalMappingCreateIS(sub_schurs->is_subs[i],&NtoS);CHKERRQ(ierr);
719       ierr = ISGlobalToLocalMappingApplyIS(NtoS,IS_GTOLM_DROP,change_primal,&change_primal_sub[i]);CHKERRQ(ierr);
720       ierr = ISLocalToGlobalMappingDestroy(&NtoS);CHKERRQ(ierr);
721     }
722     ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,change_primal,&change_primal_B);CHKERRQ(ierr);
723     ierr = ISLocalToGlobalMappingCreateIS(sub_schurs->is_Ej_all,&BtoS);CHKERRQ(ierr);
724     ierr = ISGlobalToLocalMappingApplyIS(BtoS,IS_GTOLM_DROP,change_primal_B,&change_primal_all);CHKERRQ(ierr);
725     ierr = ISLocalToGlobalMappingDestroy(&BtoS);CHKERRQ(ierr);
726     ierr = ISDestroy(&change_primal_B);CHKERRQ(ierr);
727     if (!sub_schurs->change) {
728       ierr = PetscCalloc1(sub_schurs->n_subs,&sub_schurs->change);CHKERRQ(ierr);
729     }
730     for (i=0;i<sub_schurs->n_subs;i++) {
731       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
732       if (!sub_schurs->change[i]) {
733         Mat change_sub;
734 
735         ierr = KSPCreate(PETSC_COMM_SELF,&sub_schurs->change[i]);CHKERRQ(ierr);
736         ierr = KSPSetType(sub_schurs->change[i],KSPPREONLY);CHKERRQ(ierr);
737         if (!sub_schurs->change_with_qr) {
738           ierr = MatGetSubMatrix(change,sub_schurs->is_subs[i],sub_schurs->is_subs[i],MAT_INITIAL_MATRIX,&change_sub);CHKERRQ(ierr);
739         } else {
740           Mat change_subt;
741           ierr = MatGetSubMatrix(change,sub_schurs->is_subs[i],sub_schurs->is_subs[i],MAT_INITIAL_MATRIX,&change_subt);CHKERRQ(ierr);
742           ierr = MatConvert(change_subt,MATSEQDENSE,MAT_INITIAL_MATRIX,&change_sub);CHKERRQ(ierr);
743           ierr = MatDestroy(&change_subt);CHKERRQ(ierr);
744         }
745         ierr = KSPSetOperators(sub_schurs->change[i],change_sub,change_sub);CHKERRQ(ierr);
746         ierr = MatDestroy(&change_sub);CHKERRQ(ierr);
747         ierr = KSPSetOptionsPrefix(sub_schurs->change[i],"sub_schurs_change_");CHKERRQ(ierr);
748       }
749     }
750     ierr = ISDestroy(&change_primal_all);CHKERRQ(ierr);
751   }
752 
753   /* Local matrix of all local Schur on subsets (transposed) */
754   if (!sub_schurs->S_Ej_all) {
755     ierr = MatCreate(PETSC_COMM_SELF,&sub_schurs->S_Ej_all);CHKERRQ(ierr);
756     ierr = MatSetSizes(sub_schurs->S_Ej_all,PETSC_DECIDE,PETSC_DECIDE,local_size,local_size);CHKERRQ(ierr);
757     ierr = MatSetType(sub_schurs->S_Ej_all,MATAIJ);CHKERRQ(ierr);
758     ierr = MatSeqAIJSetPreallocation(sub_schurs->S_Ej_all,0,nnz);CHKERRQ(ierr);
759   }
760 
761   /* Compute Schur complements explicitly */
762   F = NULL;
763   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 */
764     Mat         S_Ej_expl;
765     PetscScalar *work;
766     PetscInt    j,*dummy_idx;
767     PetscBool   Sdense;
768 
769     ierr = PetscMalloc2(max_subset_size,&dummy_idx,max_subset_size*max_subset_size,&work);CHKERRQ(ierr);
770     local_size = 0;
771     for (i=0;i<sub_schurs->n_subs;i++) {
772       IS  is_subset_B;
773       Mat AE_EE,AE_IE,AE_EI,S_Ej;
774 
775       /* subsets in original and boundary numbering */
776       ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,sub_schurs->is_subs[i],&is_subset_B);CHKERRQ(ierr);
777       /* EE block */
778       ierr = MatGetSubMatrix(A_BB,is_subset_B,is_subset_B,MAT_INITIAL_MATRIX,&AE_EE);CHKERRQ(ierr);
779       /* IE block */
780       ierr = MatGetSubMatrix(A_IB,is_I,is_subset_B,MAT_INITIAL_MATRIX,&AE_IE);CHKERRQ(ierr);
781       /* EI block */
782       if (sub_schurs->is_hermitian) {
783         ierr = MatCreateTranspose(AE_IE,&AE_EI);CHKERRQ(ierr);
784       } else {
785         ierr = MatGetSubMatrix(A_BI,is_subset_B,is_I,MAT_INITIAL_MATRIX,&AE_EI);CHKERRQ(ierr);
786       }
787       ierr = ISDestroy(&is_subset_B);CHKERRQ(ierr);
788       ierr = MatCreateSchurComplement(AE_II,AE_II,AE_IE,AE_EI,AE_EE,&S_Ej);CHKERRQ(ierr);
789       ierr = MatDestroy(&AE_EE);CHKERRQ(ierr);
790       ierr = MatDestroy(&AE_IE);CHKERRQ(ierr);
791       ierr = MatDestroy(&AE_EI);CHKERRQ(ierr);
792       if (AE_II == A_II) { /* we can reuse the same ksp */
793         KSP ksp;
794         ierr = MatSchurComplementGetKSP(sub_schurs->S,&ksp);CHKERRQ(ierr);
795         ierr = MatSchurComplementSetKSP(S_Ej,ksp);CHKERRQ(ierr);
796       } else { /* build new ksp object which inherits ksp and pc types from the original one */
797         KSP       origksp,schurksp;
798         PC        origpc,schurpc;
799         KSPType   ksp_type;
800         PetscInt  n_internal;
801         PetscBool ispcnone;
802 
803         ierr = MatSchurComplementGetKSP(sub_schurs->S,&origksp);CHKERRQ(ierr);
804         ierr = MatSchurComplementGetKSP(S_Ej,&schurksp);CHKERRQ(ierr);
805         ierr = KSPGetType(origksp,&ksp_type);CHKERRQ(ierr);
806         ierr = KSPSetType(schurksp,ksp_type);CHKERRQ(ierr);
807         ierr = KSPGetPC(schurksp,&schurpc);CHKERRQ(ierr);
808         ierr = KSPGetPC(origksp,&origpc);CHKERRQ(ierr);
809         ierr = PetscObjectTypeCompare((PetscObject)origpc,PCNONE,&ispcnone);CHKERRQ(ierr);
810         if (!ispcnone) {
811           PCType pc_type;
812           ierr = PCGetType(origpc,&pc_type);CHKERRQ(ierr);
813           ierr = PCSetType(schurpc,pc_type);CHKERRQ(ierr);
814         } else {
815           ierr = PCSetType(schurpc,PCLU);CHKERRQ(ierr);
816         }
817         ierr = ISGetSize(is_I,&n_internal);CHKERRQ(ierr);
818         if (n_internal) { /* UMFPACK gives error with 0 sized problems */
819           MatSolverPackage solver=NULL;
820           ierr = PCFactorGetMatSolverPackage(origpc,(const MatSolverPackage*)&solver);CHKERRQ(ierr);
821           if (solver) {
822             ierr = PCFactorSetMatSolverPackage(schurpc,solver);CHKERRQ(ierr);
823           }
824         }
825         ierr = KSPSetUp(schurksp);CHKERRQ(ierr);
826       }
827       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
828       ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&S_Ej_expl);CHKERRQ(ierr);
829       ierr = PCBDDCComputeExplicitSchur(S_Ej,sub_schurs->is_hermitian,MAT_REUSE_MATRIX,&S_Ej_expl);CHKERRQ(ierr);
830       ierr = PetscObjectTypeCompare((PetscObject)S_Ej_expl,MATSEQDENSE,&Sdense);CHKERRQ(ierr);
831       if (Sdense) {
832         for (j=0;j<subset_size;j++) {
833           dummy_idx[j]=local_size+j;
834         }
835         ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
836       } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not yet implemented for sparse matrices");
837       ierr = MatDestroy(&S_Ej);CHKERRQ(ierr);
838       ierr = MatDestroy(&S_Ej_expl);CHKERRQ(ierr);
839       local_size += subset_size;
840     }
841     ierr = PetscFree2(dummy_idx,work);CHKERRQ(ierr);
842     /* free */
843     ierr = ISDestroy(&is_I);CHKERRQ(ierr);
844     ierr = MatDestroy(&AE_II);CHKERRQ(ierr);
845     ierr = PetscFree(all_local_idx_N);CHKERRQ(ierr);
846   } else {
847     Mat         A,cs_AIB_mat = NULL,benign_AIIm1_ones_mat = NULL;
848     Vec         Dall = NULL;
849     IS          is_A_all,*is_p_r = NULL;
850     PetscScalar *work,*S_data,*schur_factor,infty = PETSC_MAX_REAL;
851     PetscInt    n,n_I,*dummy_idx,size_schur,size_active_schur,cum,cum2;
852     PetscBool   economic,solver_S,S_lower_triangular = PETSC_FALSE,factor_workaround;
853     char        solver[256];
854 
855     /* get sizes */
856     n_I = 0;
857     if (is_I_layer) {
858       ierr = ISGetLocalSize(is_I_layer,&n_I);CHKERRQ(ierr);
859     }
860     economic = PETSC_FALSE;
861     ierr = ISGetLocalSize(sub_schurs->is_I,&cum);CHKERRQ(ierr);
862     if (cum != n_I) economic = PETSC_TRUE;
863     ierr = MatGetLocalSize(sub_schurs->A,&n,NULL);CHKERRQ(ierr);
864     size_active_schur = local_size;
865 
866     /* import scaling vector */
867     if (scaling && compute_Stilda) {
868       const PetscScalar *array;
869       PetscScalar       *array2;
870       const PetscInt    *idxs;
871       PetscInt          i;
872 
873       ierr = ISGetIndices(sub_schurs->is_Ej_all,&idxs);CHKERRQ(ierr);
874       ierr = VecCreateSeq(PETSC_COMM_SELF,size_active_schur,&Dall);CHKERRQ(ierr);
875       ierr = VecGetArrayRead(scaling,&array);CHKERRQ(ierr);
876       ierr = VecGetArray(Dall,&array2);CHKERRQ(ierr);
877       for (i=0;i<size_active_schur;i++) array2[i] = array[idxs[i]];
878       ierr = VecRestoreArray(Dall,&array2);CHKERRQ(ierr);
879       ierr = VecRestoreArrayRead(scaling,&array);CHKERRQ(ierr);
880       ierr = ISRestoreIndices(sub_schurs->is_Ej_all,&idxs);CHKERRQ(ierr);
881       deluxe = PETSC_FALSE;
882     }
883 
884     /* size active schurs does not count any dirichlet or vertex dof on the interface */
885     cum = n_I+size_active_schur;
886     if (sub_schurs->is_dir) {
887       const PetscInt* idxs;
888       PetscInt        n_dir;
889 
890       ierr = ISGetLocalSize(sub_schurs->is_dir,&n_dir);CHKERRQ(ierr);
891       ierr = ISGetIndices(sub_schurs->is_dir,&idxs);CHKERRQ(ierr);
892       ierr = PetscMemcpy(all_local_idx_N+cum,idxs,n_dir*sizeof(PetscInt));CHKERRQ(ierr);
893       ierr = ISRestoreIndices(sub_schurs->is_dir,&idxs);CHKERRQ(ierr);
894       cum += n_dir;
895     }
896     factor_workaround = PETSC_FALSE;
897     /* include the primal vertices in the Schur complement */
898     if (exact_schur && sub_schurs->is_vertices && (compute_Stilda || benign_n)) {
899       PetscInt n_v;
900 
901       ierr = ISGetLocalSize(sub_schurs->is_vertices,&n_v);CHKERRQ(ierr);
902       if (n_v) {
903         const PetscInt* idxs;
904 
905         ierr = ISGetIndices(sub_schurs->is_vertices,&idxs);CHKERRQ(ierr);
906         ierr = PetscMemcpy(all_local_idx_N+cum,idxs,n_v*sizeof(PetscInt));CHKERRQ(ierr);
907         ierr = ISRestoreIndices(sub_schurs->is_vertices,&idxs);CHKERRQ(ierr);
908         cum += n_v;
909         factor_workaround = PETSC_TRUE;
910       }
911     }
912     size_schur = cum - n_I;
913     ierr = ISCreateGeneral(PETSC_COMM_SELF,cum,all_local_idx_N,PETSC_OWN_POINTER,&is_A_all);CHKERRQ(ierr);
914     /* get working mat (TODO: factorize without actually permuting)  */
915     if (cum == n) {
916       ierr = ISSetPermutation(is_A_all);CHKERRQ(ierr);
917       ierr = MatPermute(sub_schurs->A,is_A_all,is_A_all,&A);CHKERRQ(ierr);
918     } else {
919       ierr = MatGetSubMatrix(sub_schurs->A,is_A_all,is_A_all,MAT_INITIAL_MATRIX,&A);CHKERRQ(ierr);
920     }
921     ierr = MatSetOptionsPrefix(A,"sub_schurs_");CHKERRQ(ierr);
922 
923     /* if we actually change the basis for the pressures, LDL^T factors will use a lot of memory
924        n^2 instead of n^1.5 or something. This is a workaround */
925     if (benign_n) {
926       Vec                    v;
927       ISLocalToGlobalMapping N_to_reor;
928       IS                     is_p0,is_p0_p;
929       PetscScalar            *cs_AIB,*AIIm1_data;
930       PetscInt               sizeA;
931 
932       ierr = ISLocalToGlobalMappingCreateIS(is_A_all,&N_to_reor);CHKERRQ(ierr);
933       ierr = ISCreateGeneral(PETSC_COMM_SELF,benign_n,benign_p0_lidx,PETSC_COPY_VALUES,&is_p0);CHKERRQ(ierr);
934       ierr = ISGlobalToLocalMappingApplyIS(N_to_reor,IS_GTOLM_DROP,is_p0,&is_p0_p);CHKERRQ(ierr);
935       ierr = ISDestroy(&is_p0);CHKERRQ(ierr);
936       ierr = MatCreateVecs(A,&v,NULL);CHKERRQ(ierr);
937       ierr = VecGetSize(v,&sizeA);CHKERRQ(ierr);
938       ierr = MatCreateSeqDense(PETSC_COMM_SELF,sizeA,benign_n,NULL,&benign_AIIm1_ones_mat);CHKERRQ(ierr);
939       ierr = MatCreateSeqDense(PETSC_COMM_SELF,benign_n,size_schur,NULL,&cs_AIB_mat);CHKERRQ(ierr);
940       ierr = MatDenseGetArray(cs_AIB_mat,&cs_AIB);CHKERRQ(ierr);
941       ierr = MatDenseGetArray(benign_AIIm1_ones_mat,&AIIm1_data);CHKERRQ(ierr);
942       ierr = PetscMalloc1(benign_n,&is_p_r);CHKERRQ(ierr);
943       /* compute colsum of A_IB restricted to pressures */
944       for (i=0;i<benign_n;i++) {
945         Vec            benign_AIIm1_ones;
946         PetscScalar    *array;
947         const PetscInt *idxs;
948         PetscInt       j,nz;
949 
950         ierr = ISGlobalToLocalMappingApplyIS(N_to_reor,IS_GTOLM_DROP,benign_zerodiag_subs[i],&is_p_r[i]);CHKERRQ(ierr);
951         ierr = ISGetLocalSize(is_p_r[i],&nz);CHKERRQ(ierr);
952         ierr = ISGetIndices(is_p_r[i],&idxs);CHKERRQ(ierr);
953         for (j=0;j<nz;j++) AIIm1_data[idxs[j]+sizeA*i] = 1.;
954         ierr = ISRestoreIndices(is_p_r[i],&idxs);CHKERRQ(ierr);
955         ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,sizeA,AIIm1_data+sizeA*i,&benign_AIIm1_ones);CHKERRQ(ierr);
956         ierr = MatMult(A,benign_AIIm1_ones,v);CHKERRQ(ierr);
957         ierr = VecDestroy(&benign_AIIm1_ones);CHKERRQ(ierr);
958         ierr = VecGetArray(v,&array);CHKERRQ(ierr);
959         for (j=0;j<size_schur;j++) cs_AIB[i*size_schur + j] = array[j+n_I]/nz;
960         ierr = VecRestoreArray(v,&array);CHKERRQ(ierr);
961       }
962       ierr = MatDenseRestoreArray(cs_AIB_mat,&cs_AIB);CHKERRQ(ierr);
963       ierr = MatDenseRestoreArray(benign_AIIm1_ones_mat,&AIIm1_data);CHKERRQ(ierr);
964       ierr = VecDestroy(&v);CHKERRQ(ierr);
965       ierr = MatSetOption(A,MAT_KEEP_NONZERO_PATTERN,PETSC_FALSE);CHKERRQ(ierr);
966       ierr = MatZeroRowsColumnsIS(A,is_p0_p,1.0,NULL,NULL);CHKERRQ(ierr);
967       ierr = ISDestroy(&is_p0_p);CHKERRQ(ierr);
968       ierr = ISLocalToGlobalMappingDestroy(&N_to_reor);CHKERRQ(ierr);
969     }
970     ierr = MatSetOption(A,MAT_SYMMETRIC,sub_schurs->is_hermitian);CHKERRQ(ierr);
971     ierr = MatSetOption(A,MAT_HERMITIAN,sub_schurs->is_hermitian);CHKERRQ(ierr);
972     ierr = MatSetOption(A,MAT_SPD,sub_schurs->is_posdef);CHKERRQ(ierr);
973 #if defined(PETSC_HAVE_MUMPS)
974     ierr = PetscStrcpy(solver,MATSOLVERMUMPS);CHKERRQ(ierr);
975 #else
976     ierr = PetscStrcpy(solver,MATSOLVERMKL_PARDISO);CHKERRQ(ierr);
977 #endif
978     ierr = PetscOptionsGetString(NULL,"sub_schurs_","-mat_solver_package",solver,256,NULL);CHKERRQ(ierr);
979 
980     /* when using the benign subspace trick, the local Schur complements are SPD */
981     if (benign_trick) sub_schurs->is_posdef = PETSC_TRUE;
982 
983     if (n_I) { /* TODO add ordering from options */
984       IS is_schur;
985 
986       if (sub_schurs->is_hermitian && sub_schurs->is_posdef) {
987         ierr = MatGetFactor(A,solver,MAT_FACTOR_CHOLESKY,&F);CHKERRQ(ierr);
988       } else {
989         ierr = MatGetFactor(A,solver,MAT_FACTOR_LU,&F);CHKERRQ(ierr);
990       }
991       /* subsets ordered last */
992       ierr = ISCreateStride(PETSC_COMM_SELF,size_schur,n_I,1,&is_schur);CHKERRQ(ierr);
993       ierr = MatFactorSetSchurIS(F,is_schur);CHKERRQ(ierr);
994       ierr = ISDestroy(&is_schur);CHKERRQ(ierr);
995 
996       /* factorization step */
997       if (sub_schurs->is_hermitian && sub_schurs->is_posdef) {
998         ierr = MatCholeskyFactorSymbolic(F,A,NULL,NULL);CHKERRQ(ierr);
999 #if defined(PETSC_HAVE_MUMPS) /* be sure that icntl 19 is not set by command line */
1000         ierr = MatMumpsSetIcntl(F,19,2);CHKERRQ(ierr);
1001 #endif
1002         if (sub_schurs->is_posdef) {
1003           ierr = MatFactorSetSchurComplementSolverType(F,1);CHKERRQ(ierr);
1004         } else {
1005           ierr = MatFactorSetSchurComplementSolverType(F,2);CHKERRQ(ierr);
1006         }
1007         ierr = MatCholeskyFactorNumeric(F,A,NULL);CHKERRQ(ierr);
1008         S_lower_triangular = PETSC_TRUE;
1009       } else {
1010         ierr = MatLUFactorSymbolic(F,A,NULL,NULL,NULL);CHKERRQ(ierr);
1011 #if defined(PETSC_HAVE_MUMPS) /* be sure that icntl 19 is not set by command line */
1012         ierr = MatMumpsSetIcntl(F,19,3);CHKERRQ(ierr);
1013 #endif
1014         ierr = MatLUFactorNumeric(F,A,NULL);CHKERRQ(ierr);
1015       }
1016 
1017       /* get explicit Schur Complement computed during numeric factorization */
1018       ierr = MatFactorGetSchurComplement(F,&S_all);CHKERRQ(ierr);
1019       /* we can reuse the solvers if we are not using the economic version */
1020       reuse_solvers = (PetscBool)(reuse_solvers && !economic);
1021       factor_workaround = (PetscBool)(reuse_solvers && factor_workaround);
1022       solver_S = PETSC_TRUE;
1023 
1024       /* update the Schur complement with the change of basis on the pressures */
1025       if (benign_n) {
1026         PetscScalar *S_data,*cs_AIB,*AIIm1_data;
1027         Vec         v;
1028         PetscInt    sizeA;
1029 
1030         ierr = MatDenseGetArray(S_all,&S_data);CHKERRQ(ierr);
1031         ierr = MatCreateVecs(A,&v,NULL);CHKERRQ(ierr);
1032         ierr = VecGetSize(v,&sizeA);CHKERRQ(ierr);
1033 #if defined(PETSC_HAVE_MUMPS)
1034         ierr = MatMumpsSetIcntl(F,26,0);CHKERRQ(ierr);
1035 #endif
1036 #if defined(PETSC_HAVE_MKL_PARDISO)
1037         ierr = MatMkl_PardisoSetCntl(F,70,1);CHKERRQ(ierr);
1038 #endif
1039         ierr = MatDenseGetArray(cs_AIB_mat,&cs_AIB);CHKERRQ(ierr);
1040         ierr = MatDenseGetArray(benign_AIIm1_ones_mat,&AIIm1_data);CHKERRQ(ierr);
1041         for (i=0;i<benign_n;i++) {
1042           Vec            benign_AIIm1_ones;
1043           PetscScalar    *array,sum = 0.;
1044           const PetscInt *idxs;
1045           PetscInt       j,nz;
1046 
1047           ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,sizeA,AIIm1_data+sizeA*i,&benign_AIIm1_ones);CHKERRQ(ierr);
1048           ierr = VecCopy(benign_AIIm1_ones,v);CHKERRQ(ierr);
1049           ierr = MatSolve(F,v,benign_AIIm1_ones);CHKERRQ(ierr);
1050           ierr = ISGetLocalSize(is_p_r[i],&nz);CHKERRQ(ierr);
1051           ierr = ISGetIndices(is_p_r[i],&idxs);CHKERRQ(ierr);
1052           for (j=0;j<nz;j++) sum += AIIm1_data[idxs[j]+sizeA*i];
1053           sum -= 1.; /* p0 dof (eliminated) is excluded from the sum */
1054           ierr = ISRestoreIndices(is_p_r[i],&idxs);CHKERRQ(ierr);
1055           ierr = MatMult(A,benign_AIIm1_ones,v);CHKERRQ(ierr);
1056           ierr = VecGetArray(v,&array);CHKERRQ(ierr);
1057           /* perform sparse rank-1 updates on symmetric Schur (TODO: use dense linear algebra?)*/
1058           /* cs_AIB already scaled by 1./nz */;
1059           sum  = -sum;
1060           ierr = SparseRankOneUpdate(S_data,size_schur,sum,cs_AIB+i*size_schur,cs_AIB+i*size_schur);CHKERRQ(ierr);
1061           sum  = 1.;
1062           ierr = SparseRankOneUpdate(S_data,size_schur,sum,array+n_I,cs_AIB+i*size_schur);CHKERRQ(ierr);
1063           ierr = SparseRankOneUpdate(S_data,size_schur,sum,cs_AIB+i*size_schur,array+n_I);CHKERRQ(ierr);
1064           ierr = VecRestoreArray(v,&array);CHKERRQ(ierr);
1065           /* set p0 entry of AIIm1_ones to zero */
1066           ierr = ISGetIndices(is_p_r[i],&idxs);CHKERRQ(ierr);
1067           AIIm1_data[idxs[nz-1]+sizeA*i] = 0.;
1068           ierr = ISRestoreIndices(is_p_r[i],&idxs);CHKERRQ(ierr);
1069           ierr = VecDestroy(&benign_AIIm1_ones);CHKERRQ(ierr);
1070         }
1071   /* restore defaults */
1072 #if defined(PETSC_HAVE_MUMPS)
1073         ierr = MatMumpsSetIcntl(F,26,-1);CHKERRQ(ierr);
1074 #endif
1075 #if defined(PETSC_HAVE_MKL_PARDISO)
1076         ierr = MatMkl_PardisoSetCntl(F,70,0);CHKERRQ(ierr);
1077 #endif
1078         ierr = MatDenseRestoreArray(cs_AIB_mat,&cs_AIB);CHKERRQ(ierr);
1079         ierr = MatDenseRestoreArray(benign_AIIm1_ones_mat,&AIIm1_data);CHKERRQ(ierr);
1080         ierr = VecDestroy(&v);CHKERRQ(ierr);
1081         ierr = MatDenseRestoreArray(S_all,&S_data);CHKERRQ(ierr);
1082       }
1083       if (!reuse_solvers) {
1084         for (i=0;i<benign_n;i++) {
1085           ierr = ISDestroy(&is_p_r[i]);CHKERRQ(ierr);
1086         }
1087         ierr = PetscFree(is_p_r);CHKERRQ(ierr);
1088         ierr = MatDestroy(&cs_AIB_mat);CHKERRQ(ierr);
1089         ierr = MatDestroy(&benign_AIIm1_ones_mat);CHKERRQ(ierr);
1090       }
1091     } else { /* we can't use MatFactor when size_schur == size_of_the_problem */
1092       ierr = MatConvert(A,MATSEQDENSE,MAT_INITIAL_MATRIX,&S_all);CHKERRQ(ierr);
1093       reuse_solvers = PETSC_FALSE; /* TODO: why we can't reuse the solvers here? */
1094       solver_S = PETSC_FALSE;
1095     }
1096 
1097     if (reuse_solvers) {
1098       Mat                A_II,Afake;
1099       Vec                vec1_B;
1100       PCBDDCReuseSolvers msolv_ctx;
1101       PetscInt           n_R;
1102 
1103       if (sub_schurs->reuse_solver) {
1104         ierr = PCBDDCReuseSolversReset(sub_schurs->reuse_solver);CHKERRQ(ierr);
1105       } else {
1106         ierr = PetscNew(&sub_schurs->reuse_solver);CHKERRQ(ierr);
1107       }
1108       msolv_ctx = sub_schurs->reuse_solver;
1109       ierr = MatSchurComplementGetSubMatrices(sub_schurs->S,&A_II,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
1110       ierr = PetscObjectReference((PetscObject)F);CHKERRQ(ierr);
1111       msolv_ctx->F = F;
1112       ierr = MatCreateVecs(F,&msolv_ctx->sol,NULL);CHKERRQ(ierr);
1113       /* currently PETSc has no support for MatSolve(F,x,x), so cheat and let rhs and sol share the same memory */
1114       {
1115         PetscScalar *array;
1116         PetscInt    n;
1117 
1118         ierr = VecGetLocalSize(msolv_ctx->sol,&n);CHKERRQ(ierr);
1119         ierr = VecGetArray(msolv_ctx->sol,&array);CHKERRQ(ierr);
1120         ierr = VecCreateSeqWithArray(PetscObjectComm((PetscObject)msolv_ctx->sol),1,n,array,&msolv_ctx->rhs);CHKERRQ(ierr);
1121         ierr = VecRestoreArray(msolv_ctx->sol,&array);CHKERRQ(ierr);
1122       }
1123       msolv_ctx->has_vertices = factor_workaround;
1124 
1125       /* interior solver */
1126       ierr = PCCreate(PetscObjectComm((PetscObject)A_II),&msolv_ctx->interior_solver);CHKERRQ(ierr);
1127       ierr = PCSetOperators(msolv_ctx->interior_solver,A_II,A_II);CHKERRQ(ierr);
1128       ierr = PCSetType(msolv_ctx->interior_solver,PCSHELL);CHKERRQ(ierr);
1129       ierr = PCShellSetContext(msolv_ctx->interior_solver,msolv_ctx);CHKERRQ(ierr);
1130       ierr = PCShellSetApply(msolv_ctx->interior_solver,PCBDDCReuseSolvers_Interior);CHKERRQ(ierr);
1131       ierr = PCShellSetApplyTranspose(msolv_ctx->interior_solver,PCBDDCReuseSolvers_InteriorTranspose);CHKERRQ(ierr);
1132 
1133       /* correction solver */
1134       ierr = PCCreate(PetscObjectComm((PetscObject)A_II),&msolv_ctx->correction_solver);CHKERRQ(ierr);
1135       ierr = PCSetType(msolv_ctx->correction_solver,PCSHELL);CHKERRQ(ierr);
1136       ierr = PCShellSetContext(msolv_ctx->correction_solver,msolv_ctx);CHKERRQ(ierr);
1137       ierr = PCShellSetApply(msolv_ctx->correction_solver,PCBDDCReuseSolvers_Correction);CHKERRQ(ierr);
1138       ierr = PCShellSetApplyTranspose(msolv_ctx->correction_solver,PCBDDCReuseSolvers_CorrectionTranspose);CHKERRQ(ierr);
1139 
1140       /* scatter and vecs for Schur complement solver */
1141       ierr = MatCreateVecs(S_all,&msolv_ctx->sol_B,&msolv_ctx->rhs_B);CHKERRQ(ierr);
1142       ierr = MatCreateVecs(sub_schurs->S,&vec1_B,NULL);CHKERRQ(ierr);
1143       if (!factor_workaround) {
1144         ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,is_A_all,&msolv_ctx->is_B);CHKERRQ(ierr);
1145         ierr = VecScatterCreate(vec1_B,msolv_ctx->is_B,msolv_ctx->sol_B,NULL,&msolv_ctx->correction_scatter_B);CHKERRQ(ierr);
1146         ierr = PetscObjectReference((PetscObject)is_A_all);CHKERRQ(ierr);
1147         msolv_ctx->is_R = is_A_all;
1148       } else {
1149         IS              is_B_all;
1150         const PetscInt* idxs;
1151         PetscInt        dual,n_v,n;
1152 
1153         ierr = ISGetLocalSize(sub_schurs->is_vertices,&n_v);CHKERRQ(ierr);
1154         dual = size_schur - n_v;
1155         ierr = ISGetLocalSize(is_A_all,&n);CHKERRQ(ierr);
1156         ierr = ISGetIndices(is_A_all,&idxs);CHKERRQ(ierr);
1157         ierr = ISCreateGeneral(PetscObjectComm((PetscObject)is_A_all),dual,idxs+n_I,PETSC_COPY_VALUES,&is_B_all);CHKERRQ(ierr);
1158         ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,is_B_all,&msolv_ctx->is_B);CHKERRQ(ierr);
1159         ierr = ISDestroy(&is_B_all);CHKERRQ(ierr);
1160         ierr = ISCreateStride(PetscObjectComm((PetscObject)is_A_all),dual,0,1,&is_B_all);CHKERRQ(ierr);
1161         ierr = VecScatterCreate(vec1_B,msolv_ctx->is_B,msolv_ctx->sol_B,is_B_all,&msolv_ctx->correction_scatter_B);CHKERRQ(ierr);
1162         ierr = ISDestroy(&is_B_all);CHKERRQ(ierr);
1163         ierr = ISCreateGeneral(PetscObjectComm((PetscObject)is_A_all),n-n_v,idxs,PETSC_COPY_VALUES,&msolv_ctx->is_R);CHKERRQ(ierr);
1164         ierr = ISRestoreIndices(is_A_all,&idxs);CHKERRQ(ierr);
1165       }
1166       ierr = ISGetLocalSize(msolv_ctx->is_R,&n_R);CHKERRQ(ierr);
1167       ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,n_R,n_R,0,NULL,&Afake);CHKERRQ(ierr);
1168       ierr = PCSetOperators(msolv_ctx->correction_solver,Afake,Afake);CHKERRQ(ierr);
1169       ierr = MatDestroy(&Afake);CHKERRQ(ierr);
1170       ierr = VecDestroy(&vec1_B);CHKERRQ(ierr);
1171 
1172       /* communicate benign info to solver context */
1173       if (benign_n) {
1174         PetscScalar *array;
1175 
1176         msolv_ctx->benign_n = benign_n;
1177         msolv_ctx->benign_zerodiag_subs = is_p_r;
1178         ierr = PetscMalloc1(benign_n,&msolv_ctx->benign_save_vals);CHKERRQ(ierr);
1179         msolv_ctx->benign_csAIB = cs_AIB_mat;
1180         ierr = MatCreateVecs(cs_AIB_mat,NULL,&msolv_ctx->benign_corr_work);CHKERRQ(ierr);
1181         ierr = VecGetArray(msolv_ctx->benign_corr_work,&array);CHKERRQ(ierr);
1182         ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,size_schur,array,&msolv_ctx->benign_dummy_schur_vec);CHKERRQ(ierr);
1183         ierr = VecRestoreArray(msolv_ctx->benign_corr_work,&array);CHKERRQ(ierr);
1184         msolv_ctx->benign_AIIm1ones = benign_AIIm1_ones_mat;
1185       }
1186     }
1187     ierr = MatDestroy(&A);CHKERRQ(ierr);
1188     ierr = ISDestroy(&is_A_all);CHKERRQ(ierr);
1189 
1190     /* Work arrays */
1191     ierr = PetscMalloc2(max_subset_size,&dummy_idx,max_subset_size*max_subset_size,&work);CHKERRQ(ierr);
1192 
1193     /* matrices for adaptive selection */
1194     if (compute_Stilda) {
1195       if (!sub_schurs->sum_S_Ej_tilda_all) {
1196         ierr = MatCreate(PETSC_COMM_SELF,&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
1197         ierr = MatSetSizes(sub_schurs->sum_S_Ej_tilda_all,PETSC_DECIDE,PETSC_DECIDE,size_active_schur,size_active_schur);CHKERRQ(ierr);
1198         ierr = MatSetType(sub_schurs->sum_S_Ej_tilda_all,MATAIJ);CHKERRQ(ierr);
1199         ierr = MatSeqAIJSetPreallocation(sub_schurs->sum_S_Ej_tilda_all,0,nnz);CHKERRQ(ierr);
1200       }
1201       if (!sub_schurs->sum_S_Ej_inv_all && deluxe) {
1202         ierr = MatCreate(PETSC_COMM_SELF,&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
1203         ierr = MatSetSizes(sub_schurs->sum_S_Ej_inv_all,PETSC_DECIDE,PETSC_DECIDE,size_active_schur,size_active_schur);CHKERRQ(ierr);
1204         ierr = MatSetType(sub_schurs->sum_S_Ej_inv_all,MATAIJ);CHKERRQ(ierr);
1205         ierr = MatSeqAIJSetPreallocation(sub_schurs->sum_S_Ej_inv_all,0,nnz);CHKERRQ(ierr);
1206       }
1207     }
1208 
1209     /* S_Ej_all */
1210     cum = cum2 = 0;
1211     ierr = MatDenseGetArray(S_all,&S_data);CHKERRQ(ierr);
1212     for (i=0;i<sub_schurs->n_subs;i++) {
1213       PetscInt j;
1214 
1215       /* get S_E */
1216       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
1217       if (S_lower_triangular) { /* I need to expand the upper triangular data (column oriented) */
1218         PetscInt k;
1219         for (k=0;k<subset_size;k++) {
1220           for (j=k;j<subset_size;j++) {
1221             work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
1222             work[j*subset_size+k] = PetscConj(S_data[cum2+k*size_schur+j]);
1223           }
1224         }
1225       } else { /* just copy to workspace */
1226         PetscInt k;
1227         for (k=0;k<subset_size;k++) {
1228           for (j=0;j<subset_size;j++) {
1229             work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
1230           }
1231         }
1232       }
1233       /* insert S_E values */
1234       for (j=0;j<subset_size;j++) dummy_idx[j] = cum+j;
1235       if (change) {
1236         Mat            change_sub,SEj,T;
1237         const PetscInt *idxs;
1238         PetscInt       nz;
1239 
1240         /* change basis */
1241         ierr = KSPGetOperators(sub_schurs->change[i],&change_sub,NULL);CHKERRQ(ierr);
1242         ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&SEj);CHKERRQ(ierr);
1243         if (!sub_schurs->change_with_qr) { /* currently there's no support for PtAP with P SeqAIJ */
1244           Mat T2;
1245           ierr = MatTransposeMatMult(change_sub,SEj,MAT_INITIAL_MATRIX,1.0,&T2);CHKERRQ(ierr);
1246           ierr = MatMatMult(T2,change_sub,MAT_INITIAL_MATRIX,1.0,&T);CHKERRQ(ierr);
1247           ierr = MatConvert(T,MATSEQDENSE,MAT_INPLACE_MATRIX,&T);CHKERRQ(ierr);
1248           ierr = MatDestroy(&T2);CHKERRQ(ierr);
1249         } else {
1250           ierr = MatPtAP(SEj,change_sub,MAT_INITIAL_MATRIX,1.0,&T);CHKERRQ(ierr);
1251         }
1252         ierr = MatCopy(T,SEj,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
1253         ierr = MatDestroy(&T);CHKERRQ(ierr);
1254         /* currently there's no support for ZeroRowsColumns with A dense */
1255         ierr = ISGetIndices(change_primal_sub[i],&idxs);CHKERRQ(ierr);
1256         ierr = ISGetLocalSize(change_primal_sub[i],&nz);CHKERRQ(ierr);
1257         for (j=0;j<nz;j++) {
1258           PetscInt k;
1259 
1260           ierr = PetscMemzero(work+subset_size*idxs[j],subset_size*sizeof(PetscScalar));CHKERRQ(ierr);
1261           for (k=0;k<subset_size;k++) work[idxs[j]+k*subset_size] = 0.;
1262           work[idxs[j]*(subset_size+1)] = 1.0;
1263         }
1264         ierr = ISRestoreIndices(change_primal_sub[i],&idxs);CHKERRQ(ierr);
1265         ierr = MatDestroy(&SEj);CHKERRQ(ierr);
1266       }
1267       if (deluxe) {
1268         ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
1269         /* if adaptivity is requested, invert S_E blocks */
1270         if (compute_Stilda) {
1271           ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr);
1272           ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
1273           if (!use_getr) { /* TODO add sytrf/i for symmetric non hermitian */
1274             PetscInt k;
1275             PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,work,&B_N,&B_ierr));
1276             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
1277             PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,work,&B_N,&B_ierr));
1278             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
1279             for (k=0;k<subset_size;k++) {
1280               for (j=k;j<subset_size;j++) {
1281                 work[j*subset_size+k] = work[k*subset_size+j];
1282               }
1283             }
1284           } else {
1285             PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,work,&B_N,pivots,&B_ierr));
1286             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
1287             PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,work,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
1288             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
1289           }
1290           ierr = PetscFPTrapPop();CHKERRQ(ierr);
1291           ierr = MatSetValues(sub_schurs->sum_S_Ej_inv_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
1292         }
1293       } else if (compute_Stilda) { /* not using deluxe */
1294         Mat         SEj;
1295         Vec         D;
1296         PetscScalar *array;
1297 
1298         ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&SEj);CHKERRQ(ierr);
1299         ierr = VecGetArray(Dall,&array);CHKERRQ(ierr);
1300         ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,subset_size,array+cum,&D);CHKERRQ(ierr);
1301         ierr = VecRestoreArray(Dall,&array);CHKERRQ(ierr);
1302         ierr = MatDiagonalScale(SEj,D,D);CHKERRQ(ierr);
1303         ierr = MatDestroy(&SEj);CHKERRQ(ierr);
1304         ierr = VecDestroy(&D);CHKERRQ(ierr);
1305         ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
1306       }
1307       cum += subset_size;
1308       cum2 += subset_size*(size_schur + 1);
1309     }
1310     ierr = MatDenseRestoreArray(S_all,&S_data);CHKERRQ(ierr);
1311 
1312     if (solver_S) {
1313       ierr = MatFactorRestoreSchurComplement(F,&S_all);CHKERRQ(ierr);
1314     }
1315 
1316     schur_factor = NULL;
1317     if (compute_Stilda && size_active_schur) {
1318 
1319       if (sub_schurs->n_subs == 1 && size_schur == size_active_schur && deluxe) { /* we already computed the inverse */
1320         PetscInt j;
1321         for (j=0;j<size_schur;j++) dummy_idx[j] = j;
1322         ierr = MatSetValues(sub_schurs->sum_S_Ej_tilda_all,size_schur,dummy_idx,size_schur,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
1323       } else {
1324         Mat S_all_inv=NULL;
1325         if (solver_S) { /* use MatFactor calls to invert S */
1326             /* let MatFactor factor the Schur complement */
1327           ierr = MatFactorFactorizeSchurComplement(F);CHKERRQ(ierr);
1328           /* for adaptive selection we need S^-1; for solver reusage we need S_\Delta\Delta^-1.
1329              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 */
1330           if (factor_workaround) {
1331             PetscScalar *data;
1332             PetscInt nd = 0;
1333 
1334             if (!sub_schurs->is_posdef) {
1335               SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor update not yet implemented for non SPD matrices");
1336             }
1337             ierr = MatFactorGetSchurComplement(F,&S_all_inv);CHKERRQ(ierr);
1338             ierr = MatDenseGetArray(S_all_inv,&data);CHKERRQ(ierr);
1339             if (sub_schurs->is_dir) { /* dirichlet dofs could have different scalings */
1340               ierr = ISGetLocalSize(sub_schurs->is_dir,&nd);CHKERRQ(ierr);
1341             }
1342             ierr = PetscMalloc1((size_active_schur*(size_active_schur+1))/2+nd,&schur_factor);CHKERRQ(ierr);
1343             cum = 0;
1344             for (i=0;i<size_active_schur;i++) {
1345               ierr = PetscMemcpy(schur_factor+cum,data+i*(size_schur+1),(size_active_schur-i)*sizeof(PetscScalar));CHKERRQ(ierr);
1346               cum += size_active_schur-i;
1347             }
1348             for (i=0;i<nd;i++) schur_factor[cum+i] = data[(i+size_active_schur)*(size_schur+1)];
1349             /* invert without calling MatFactorInvertSchurComplement, since we are hacking */
1350             ierr = PetscBLASIntCast(size_schur,&B_N);CHKERRQ(ierr);
1351             ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
1352             PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,data,&B_N,&B_ierr));
1353             ierr = PetscFPTrapPop();CHKERRQ(ierr);
1354             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
1355             ierr = MatDenseRestoreArray(S_all_inv,&data);CHKERRQ(ierr);
1356           } else {
1357             ierr = MatFactorInvertSchurComplement(F);CHKERRQ(ierr);
1358             ierr = MatFactorGetSchurComplement(F,&S_all_inv);CHKERRQ(ierr);
1359           }
1360         } else { /* we need to invert explicitly since we are not using MatFactor for S */
1361           ierr = PetscObjectReference((PetscObject)S_all);CHKERRQ(ierr);
1362           S_all_inv = S_all;
1363           ierr = MatDenseGetArray(S_all_inv,&S_data);CHKERRQ(ierr);
1364           ierr = PetscBLASIntCast(size_schur,&B_N);CHKERRQ(ierr);
1365           ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
1366           if (!use_getr) { /* TODO add sytrf/i for symmetric non hermitian */
1367             PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,S_data,&B_N,&B_ierr));
1368             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
1369             PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,S_data,&B_N,&B_ierr));
1370             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
1371           } else {
1372             PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,S_data,&B_N,pivots,&B_ierr));
1373             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
1374             PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,S_data,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
1375             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
1376           }
1377           ierr = PetscFPTrapPop();CHKERRQ(ierr);
1378           ierr = MatDenseRestoreArray(S_all_inv,&S_data);CHKERRQ(ierr);
1379         }
1380         /* S_Ej_tilda_all */
1381         cum = cum2 = 0;
1382         ierr = MatDenseGetArray(S_all_inv,&S_data);CHKERRQ(ierr);
1383         for (i=0;i<sub_schurs->n_subs;i++) {
1384           PetscInt j;
1385 
1386           ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
1387           /* get (St^-1)_E */
1388           /* Unless we are changing the variables, I don't need to expand to upper triangular since St^-1
1389              will be properly accessed later during adaptive selection */
1390           if (S_lower_triangular) {
1391             PetscInt k;
1392             if (change) {
1393               for (k=0;k<subset_size;k++) {
1394                 for (j=k;j<subset_size;j++) {
1395                   work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
1396                   work[j*subset_size+k] = work[k*subset_size+j];
1397                 }
1398               }
1399             } else {
1400               for (k=0;k<subset_size;k++) {
1401                 for (j=k;j<subset_size;j++) {
1402                   work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
1403                 }
1404               }
1405             }
1406           } else {
1407             PetscInt k;
1408             for (k=0;k<subset_size;k++) {
1409               for (j=0;j<subset_size;j++) {
1410                 work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
1411               }
1412             }
1413           }
1414           if (change) {
1415             Mat change_sub,SEj,T;
1416 
1417             /* change basis */
1418             ierr = KSPGetOperators(sub_schurs->change[i],&change_sub,NULL);CHKERRQ(ierr);
1419             ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&SEj);CHKERRQ(ierr);
1420             if (!sub_schurs->change_with_qr) { /* currently there's no support for PtAP with P SeqAIJ */
1421               Mat T2;
1422               ierr = MatTransposeMatMult(change_sub,SEj,MAT_INITIAL_MATRIX,1.0,&T2);CHKERRQ(ierr);
1423               ierr = MatMatMult(T2,change_sub,MAT_INITIAL_MATRIX,1.0,&T);CHKERRQ(ierr);
1424               ierr = MatDestroy(&T2);CHKERRQ(ierr);
1425               ierr = MatConvert(T,MATSEQDENSE,MAT_INPLACE_MATRIX,&T);CHKERRQ(ierr);
1426             } else {
1427               ierr = MatPtAP(SEj,change_sub,MAT_INITIAL_MATRIX,1.0,&T);CHKERRQ(ierr);
1428             }
1429             ierr = MatCopy(T,SEj,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
1430             ierr = MatDestroy(&T);CHKERRQ(ierr);
1431             /* set diagonal entry to a very large value to pick the basis we are eliminating as the first eigenvectors with adaptive selection */
1432             {
1433               const PetscInt *idxs;
1434               PetscInt n,j;
1435               ierr = ISGetLocalSize(change_primal_sub[i],&n);CHKERRQ(ierr);
1436               ierr = ISGetIndices(change_primal_sub[i],&idxs);CHKERRQ(ierr);
1437               for (j=0;j<n;j++) {
1438                 PetscInt k;
1439                 for (k=0;k<subset_size;k++) {
1440                   work[k+idxs[j]*subset_size] = work[idxs[j]+k*subset_size] = 0.;
1441                 }
1442                 work[idxs[j] + idxs[j]*subset_size] = 1./PETSC_SMALL;
1443               }
1444               ierr = ISRestoreIndices(change_primal_sub[i],&idxs);CHKERRQ(ierr);
1445             }
1446             ierr = MatDestroy(&SEj);CHKERRQ(ierr);
1447           }
1448           for (j=0;j<subset_size;j++) dummy_idx[j] = cum+j;
1449           ierr = MatSetValues(sub_schurs->sum_S_Ej_tilda_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
1450           cum += subset_size;
1451           cum2 += subset_size*(size_schur + 1);
1452         }
1453         ierr = MatDenseRestoreArray(S_all_inv,&S_data);CHKERRQ(ierr);
1454         if (solver_S) {
1455           ierr = MatFactorRestoreSchurComplement(F,&S_all_inv);CHKERRQ(ierr);
1456         }
1457         ierr = MatDestroy(&S_all_inv);CHKERRQ(ierr);
1458       }
1459 
1460       /* move back factors to Schur data into F */
1461       if (factor_workaround) {
1462         Mat S_tmp;
1463         PetscInt nv,nd=0;
1464 
1465         if (!solver_S) {
1466           SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"This should not happen");
1467         }
1468         ierr = ISGetLocalSize(sub_schurs->is_vertices,&nv);CHKERRQ(ierr);
1469         ierr = MatFactorGetSchurComplement(F,&S_tmp);CHKERRQ(ierr);
1470         if (sub_schurs->is_posdef) {
1471           PetscScalar *data;
1472 
1473           ierr = MatZeroEntries(S_tmp);CHKERRQ(ierr);
1474           ierr = MatDenseGetArray(S_tmp,&data);CHKERRQ(ierr);
1475 
1476           if (S_lower_triangular) {
1477             cum = 0;
1478             for (i=0;i<size_active_schur;i++) {
1479               ierr = PetscMemcpy(data+i*(size_schur+1),schur_factor+cum,(size_active_schur-i)*sizeof(PetscScalar));CHKERRQ(ierr);
1480               cum += size_active_schur-i;
1481 	    }
1482           } else {
1483             ierr = PetscMemcpy(data,schur_factor,size_schur*size_schur*sizeof(PetscScalar));CHKERRQ(ierr);
1484           }
1485           if (sub_schurs->is_dir) {
1486             ierr = ISGetLocalSize(sub_schurs->is_dir,&nd);CHKERRQ(ierr);
1487 	    for (i=0;i<nd;i++) {
1488 	      data[(i+size_active_schur)*(size_schur+1)] = schur_factor[cum+i];
1489 	    }
1490           }
1491           /* workaround: since I cannot modify the matrices used inside the solvers for the forward and backward substitutions,
1492              set the diagonal entry of the Schur factor to a very large value */
1493           for (i=size_active_schur+nd;i<size_schur;i++) {
1494             data[i*(size_schur+1)] = infty;
1495           }
1496           ierr = MatDenseRestoreArray(S_tmp,&data);CHKERRQ(ierr);
1497         } else {
1498           SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor update not yet implemented for non SPD matrices");
1499         }
1500         ierr = MatFactorRestoreSchurComplement(F,&S_tmp);CHKERRQ(ierr);
1501       }
1502     } else if (factor_workaround) { /* we need to eliminate any unneeded coupling */
1503       PetscScalar *data;
1504       PetscInt    nd = 0;
1505 
1506       if (sub_schurs->is_dir) { /* dirichlet dofs could have different scalings */
1507         ierr = ISGetLocalSize(sub_schurs->is_dir,&nd);CHKERRQ(ierr);
1508       }
1509       ierr = MatFactorGetSchurComplement(F,&S_all);CHKERRQ(ierr);
1510       ierr = MatDenseGetArray(S_all,&data);CHKERRQ(ierr);
1511       for (i=0;i<size_active_schur;i++) {
1512         ierr = PetscMemzero(data+i*size_schur+size_active_schur,(size_schur-size_active_schur)*sizeof(PetscScalar));CHKERRQ(ierr);
1513       }
1514       for (i=size_active_schur+nd;i<size_schur;i++) {
1515         ierr = PetscMemzero(data+i*size_schur+size_active_schur,(size_schur-size_active_schur)*sizeof(PetscScalar));CHKERRQ(ierr);
1516         data[i*(size_schur+1)] = infty;
1517       }
1518       ierr = MatDenseRestoreArray(S_all,&data);CHKERRQ(ierr);
1519       ierr = MatFactorRestoreSchurComplement(F,&S_all);CHKERRQ(ierr);
1520     }
1521     ierr = PetscFree2(dummy_idx,work);CHKERRQ(ierr);
1522     ierr = PetscFree(schur_factor);CHKERRQ(ierr);
1523     ierr = VecDestroy(&Dall);CHKERRQ(ierr);
1524   }
1525   if (change_primal_sub) {
1526     for (i=0;i<sub_schurs->n_subs;i++) {
1527       ierr = ISDestroy(&change_primal_sub[i]);CHKERRQ(ierr);
1528     }
1529   }
1530   ierr = PetscFree(change_primal_sub);CHKERRQ(ierr);
1531   ierr = PetscFree(nnz);CHKERRQ(ierr);
1532   ierr = MatDestroy(&F);CHKERRQ(ierr);
1533   ierr = ISDestroy(&is_I_layer);CHKERRQ(ierr);
1534   ierr = MatDestroy(&S_all);CHKERRQ(ierr);
1535   ierr = MatDestroy(&A_BB);CHKERRQ(ierr);
1536   ierr = MatDestroy(&A_IB);CHKERRQ(ierr);
1537   ierr = MatDestroy(&A_BI);CHKERRQ(ierr);
1538   ierr = MatAssemblyBegin(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1539   ierr = MatAssemblyEnd(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1540   if (compute_Stilda) {
1541     ierr = MatAssemblyBegin(sub_schurs->sum_S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1542     ierr = MatAssemblyEnd(sub_schurs->sum_S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1543     if (deluxe) {
1544       ierr = MatAssemblyBegin(sub_schurs->sum_S_Ej_inv_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1545       ierr = MatAssemblyEnd(sub_schurs->sum_S_Ej_inv_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1546     }
1547   }
1548 
1549   /* Global matrix of all assembled Schur on subsets */
1550   ierr = MatISSetLocalMat(work_mat,sub_schurs->S_Ej_all);CHKERRQ(ierr);
1551   ierr = MatISSetMPIXAIJPreallocation_Private(work_mat,global_schur_subsets,PETSC_TRUE);CHKERRQ(ierr);
1552   ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
1553 
1554   /* Get local part of (\sum_j S_Ej) */
1555   if (!sub_schurs->sum_S_Ej_all) {
1556     ierr = MatGetSubMatrices(global_schur_subsets,1,&all_subsets_n,&all_subsets_n,MAT_INITIAL_MATRIX,&submats);CHKERRQ(ierr);
1557     sub_schurs->sum_S_Ej_all = submats[0];
1558   } else {
1559     ierr = PetscMalloc1(1,&submats);CHKERRQ(ierr);
1560     submats[0] = sub_schurs->sum_S_Ej_all;
1561     ierr = MatGetSubMatrices(global_schur_subsets,1,&all_subsets_n,&all_subsets_n,MAT_REUSE_MATRIX,&submats);CHKERRQ(ierr);
1562   }
1563 
1564 
1565   /* Compute explicitly (\sum_j S_Ej)^-1 (faster scaling during PCApply, needs extra work when doing setup) */
1566 #if 0
1567   if (faster_deluxe) {
1568     Mat         tmpmat;
1569     PetscScalar *array;
1570     PetscInt    cum;
1571 
1572     ierr = MatSeqAIJGetArray(sub_schurs->sum_S_Ej_all,&array);CHKERRQ(ierr);
1573     cum = 0;
1574     for (i=0;i<sub_schurs->n_subs;i++) {
1575       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
1576       ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr);
1577       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
1578       if (!use_getr) {
1579         PetscInt j,k;
1580 
1581         PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,array+cum,&B_N,&B_ierr));
1582         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
1583         PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,array+cum,&B_N,&B_ierr));
1584         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
1585         for (j=0;j<B_N;j++) {
1586           for (k=j+1;k<B_N;k++) {
1587             array[k*B_N+j+cum] = array[j*B_N+k+cum];
1588           }
1589         }
1590       } else {
1591         PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,array+cum,&B_N,pivots,&B_ierr));
1592         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
1593         PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,array+cum,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
1594         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
1595       }
1596       ierr = PetscFPTrapPop();CHKERRQ(ierr);
1597       cum += subset_size*subset_size;
1598     }
1599     ierr = MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_all,&array);CHKERRQ(ierr);
1600     ierr = MatMatMult(sub_schurs->S_Ej_all,sub_schurs->sum_S_Ej_all,MAT_INITIAL_MATRIX,1.0,&tmpmat);CHKERRQ(ierr);
1601     ierr = MatDestroy(&sub_schurs->S_Ej_all);CHKERRQ(ierr);
1602     ierr = MatDestroy(&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
1603     sub_schurs->S_Ej_all = tmpmat;
1604   }
1605 #endif
1606   /* Get local part of (\sum_j S^-1_Ej) (\sum_j St^-1_Ej) */
1607   if (compute_Stilda) {
1608     ierr = MatISSetLocalMat(work_mat,sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
1609     ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
1610     submats[0] = sub_schurs->sum_S_Ej_tilda_all;
1611     ierr = MatGetSubMatrices(global_schur_subsets,1,&all_subsets_n,&all_subsets_n,MAT_REUSE_MATRIX,&submats);CHKERRQ(ierr);
1612     if (deluxe) {
1613       ierr = MatISSetLocalMat(work_mat,sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
1614       ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
1615       submats[0] = sub_schurs->sum_S_Ej_inv_all;
1616       ierr = MatGetSubMatrices(global_schur_subsets,1,&all_subsets_n,&all_subsets_n,MAT_REUSE_MATRIX,&submats);CHKERRQ(ierr);
1617     } else {
1618       PetscScalar *array;
1619       PetscInt    cum;
1620 
1621       ierr = MatSeqAIJGetArray(sub_schurs->sum_S_Ej_tilda_all,&array);CHKERRQ(ierr);
1622       cum = 0;
1623       for (i=0;i<sub_schurs->n_subs;i++) {
1624         ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
1625         ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr);
1626         ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
1627         PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,array+cum,&B_N,&B_ierr));
1628         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
1629         PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,array+cum,&B_N,&B_ierr));
1630         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
1631         ierr = PetscFPTrapPop();CHKERRQ(ierr);
1632         cum += subset_size*subset_size;
1633       }
1634       ierr = MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_tilda_all,&array);CHKERRQ(ierr);
1635       ierr = PetscObjectReference((PetscObject)sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
1636       ierr = MatDestroy(&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
1637       sub_schurs->sum_S_Ej_inv_all = sub_schurs->sum_S_Ej_all;
1638     }
1639   }
1640 
1641   /* free workspace */
1642   ierr = PetscFree(submats);CHKERRQ(ierr);
1643   ierr = PetscFree2(Bwork,pivots);CHKERRQ(ierr);
1644   ierr = MatDestroy(&global_schur_subsets);CHKERRQ(ierr);
1645   ierr = MatDestroy(&work_mat);CHKERRQ(ierr);
1646   ierr = ISDestroy(&all_subsets_n);CHKERRQ(ierr);
1647   ierr = PetscCommDestroy(&comm_n);CHKERRQ(ierr);
1648   PetscFunctionReturn(0);
1649 }
1650 
1651 #undef __FUNCT__
1652 #define __FUNCT__ "PCBDDCSubSchursInit"
1653 PetscErrorCode PCBDDCSubSchursInit(PCBDDCSubSchurs sub_schurs, IS is_I, IS is_B, PCBDDCGraph graph, ISLocalToGlobalMapping BtoNmap)
1654 {
1655   IS              *faces,*edges,*all_cc,vertices;
1656   PetscInt        i,n_faces,n_edges,n_all_cc;
1657   PetscBool       is_sorted;
1658   PetscErrorCode  ierr;
1659 
1660   PetscFunctionBegin;
1661   ierr = ISSorted(is_I,&is_sorted);CHKERRQ(ierr);
1662   if (!is_sorted) SETERRQ(PetscObjectComm((PetscObject)is_I),PETSC_ERR_PLIB,"IS for I dofs should be shorted");
1663   ierr = ISSorted(is_B,&is_sorted);CHKERRQ(ierr);
1664   if (!is_sorted) SETERRQ(PetscObjectComm((PetscObject)is_B),PETSC_ERR_PLIB,"IS for B dofs should be shorted");
1665 
1666   /* reset any previous data */
1667   ierr = PCBDDCSubSchursReset(sub_schurs);CHKERRQ(ierr);
1668 
1669   /* get index sets for faces and edges (already sorted by global ordering) */
1670   ierr = PCBDDCGraphGetCandidatesIS(graph,&n_faces,&faces,&n_edges,&edges,&vertices);CHKERRQ(ierr);
1671   n_all_cc = n_faces+n_edges;
1672   ierr = PetscBTCreate(n_all_cc,&sub_schurs->is_edge);CHKERRQ(ierr);
1673   ierr = PetscMalloc1(n_all_cc,&all_cc);CHKERRQ(ierr);
1674   for (i=0;i<n_faces;i++) {
1675     all_cc[i] = faces[i];
1676   }
1677   for (i=0;i<n_edges;i++) {
1678     all_cc[n_faces+i] = edges[i];
1679     ierr = PetscBTSet(sub_schurs->is_edge,n_faces+i);CHKERRQ(ierr);
1680   }
1681   ierr = PetscFree(faces);CHKERRQ(ierr);
1682   ierr = PetscFree(edges);CHKERRQ(ierr);
1683   sub_schurs->is_dir = NULL;
1684   ierr = PCBDDCGraphGetDirichletDofsB(graph,&sub_schurs->is_dir);CHKERRQ(ierr);
1685 
1686   /* Determine if MatFactor can be used */
1687   sub_schurs->schur_explicit = PETSC_FALSE;
1688 #if defined(PETSC_HAVE_MUMPS)
1689   sub_schurs->schur_explicit = PETSC_TRUE;
1690 #endif
1691 #if defined(PETSC_HAVE_MKL_PARDISO)
1692   sub_schurs->schur_explicit = PETSC_TRUE;
1693 #endif
1694 
1695   ierr = PetscObjectReference((PetscObject)is_I);CHKERRQ(ierr);
1696   sub_schurs->is_I = is_I;
1697   ierr = PetscObjectReference((PetscObject)is_B);CHKERRQ(ierr);
1698   sub_schurs->is_B = is_B;
1699   ierr = PetscObjectReference((PetscObject)graph->l2gmap);CHKERRQ(ierr);
1700   sub_schurs->l2gmap = graph->l2gmap;
1701   ierr = PetscObjectReference((PetscObject)BtoNmap);CHKERRQ(ierr);
1702   sub_schurs->BtoNmap = BtoNmap;
1703   sub_schurs->n_subs = n_all_cc;
1704   sub_schurs->is_subs = all_cc;
1705   if (!sub_schurs->schur_explicit) { /* sort by local ordering if we are not using MatFactor */
1706     for (i=0;i<sub_schurs->n_subs;i++) {
1707       ierr = ISSort(sub_schurs->is_subs[i]);CHKERRQ(ierr);
1708     }
1709   }
1710   sub_schurs->is_vertices = vertices;
1711   sub_schurs->S_Ej_all = NULL;
1712   sub_schurs->sum_S_Ej_all = NULL;
1713   sub_schurs->sum_S_Ej_inv_all = NULL;
1714   sub_schurs->sum_S_Ej_tilda_all = NULL;
1715   sub_schurs->is_Ej_all = NULL;
1716   PetscFunctionReturn(0);
1717 }
1718 
1719 #undef __FUNCT__
1720 #define __FUNCT__ "PCBDDCSubSchursCreate"
1721 PetscErrorCode PCBDDCSubSchursCreate(PCBDDCSubSchurs *sub_schurs)
1722 {
1723   PCBDDCSubSchurs schurs_ctx;
1724   PetscErrorCode  ierr;
1725 
1726   PetscFunctionBegin;
1727   ierr = PetscNew(&schurs_ctx);CHKERRQ(ierr);
1728   schurs_ctx->n_subs = 0;
1729   *sub_schurs = schurs_ctx;
1730   PetscFunctionReturn(0);
1731 }
1732 
1733 #undef __FUNCT__
1734 #define __FUNCT__ "PCBDDCSubSchursReset"
1735 PetscErrorCode PCBDDCSubSchursReset(PCBDDCSubSchurs sub_schurs)
1736 {
1737   PetscInt       i;
1738   PetscErrorCode ierr;
1739 
1740   PetscFunctionBegin;
1741   ierr = MatDestroy(&sub_schurs->A);CHKERRQ(ierr);
1742   ierr = MatDestroy(&sub_schurs->S);CHKERRQ(ierr);
1743   ierr = ISDestroy(&sub_schurs->is_I);CHKERRQ(ierr);
1744   ierr = ISDestroy(&sub_schurs->is_B);CHKERRQ(ierr);
1745   ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->l2gmap);CHKERRQ(ierr);
1746   ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->BtoNmap);CHKERRQ(ierr);
1747   ierr = MatDestroy(&sub_schurs->S_Ej_all);CHKERRQ(ierr);
1748   ierr = MatDestroy(&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
1749   ierr = MatDestroy(&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
1750   ierr = MatDestroy(&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
1751   ierr = ISDestroy(&sub_schurs->is_Ej_all);CHKERRQ(ierr);
1752   ierr = ISDestroy(&sub_schurs->is_vertices);CHKERRQ(ierr);
1753   ierr = ISDestroy(&sub_schurs->is_dir);CHKERRQ(ierr);
1754   ierr = PetscBTDestroy(&sub_schurs->is_edge);CHKERRQ(ierr);
1755   for (i=0;i<sub_schurs->n_subs;i++) {
1756     ierr = ISDestroy(&sub_schurs->is_subs[i]);CHKERRQ(ierr);
1757   }
1758   if (sub_schurs->n_subs) {
1759     ierr = PetscFree(sub_schurs->is_subs);CHKERRQ(ierr);
1760   }
1761   if (sub_schurs->reuse_solver) {
1762     ierr = PCBDDCReuseSolversReset(sub_schurs->reuse_solver);CHKERRQ(ierr);
1763   }
1764   ierr = PetscFree(sub_schurs->reuse_solver);CHKERRQ(ierr);
1765   if (sub_schurs->change) {
1766     for (i=0;i<sub_schurs->n_subs;i++) {
1767       ierr = KSPDestroy(&sub_schurs->change[i]);CHKERRQ(ierr);
1768     }
1769   }
1770   ierr = PetscFree(sub_schurs->change);CHKERRQ(ierr);
1771   sub_schurs->n_subs = 0;
1772   PetscFunctionReturn(0);
1773 }
1774 
1775 #undef __FUNCT__
1776 #define __FUNCT__ "PCBDDCAdjGetNextLayer_Private"
1777 PETSC_STATIC_INLINE PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt* queue_tip,PetscInt n_prev,PetscBT touched,PetscInt* xadj,PetscInt* adjncy,PetscInt* n_added)
1778 {
1779   PetscInt       i,j,n;
1780   PetscErrorCode ierr;
1781 
1782   PetscFunctionBegin;
1783   n = 0;
1784   for (i=-n_prev;i<0;i++) {
1785     PetscInt start_dof = queue_tip[i];
1786     for (j=xadj[start_dof];j<xadj[start_dof+1];j++) {
1787       PetscInt dof = adjncy[j];
1788       if (!PetscBTLookup(touched,dof)) {
1789         ierr = PetscBTSet(touched,dof);CHKERRQ(ierr);
1790         queue_tip[n] = dof;
1791         n++;
1792       }
1793     }
1794   }
1795   *n_added = n;
1796   PetscFunctionReturn(0);
1797 }
1798