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