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