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