xref: /petsc/src/ksp/pc/impls/bddc/bddcschurs.c (revision df4d28bf1956980d02873ddcd68bc64d39668ce4)
1 #include <../src/ksp/pc/impls/bddc/bddc.h>
2 #include <../src/ksp/pc/impls/bddc/bddcprivate.h>
3 #include <petscblaslapack.h>
4 
5 PETSC_STATIC_INLINE PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt*,PetscInt,PetscBT,PetscInt*,PetscInt*,PetscInt*);
6 static PetscErrorCode PCBDDCComputeExplicitSchur(Mat,PetscBool,MatReuse,Mat*);
7 static PetscErrorCode PCBDDCReuseSolvers_Interior(PC,Vec,Vec);
8 static PetscErrorCode PCBDDCReuseSolvers_Correction(PC,Vec,Vec);
9 
10 
11 #undef __FUNCT__
12 #define __FUNCT__ "SparseRankOneUpdate"
13 PETSC_STATIC_INLINE PetscErrorCode SparseRankOneUpdate(PetscScalar *S,PetscInt n,PetscScalar a,PetscScalar *v1,PetscScalar *v2)
14 {
15   PetscInt i;
16 
17   PetscFunctionBegin;
18   for (i=0;i<n;i++) {
19     if (PetscUnlikely(PetscAbsScalar(v2[i]) > PETSC_SMALL)) {
20       PetscScalar  v2v = a*v2[i];
21       PetscBLASInt B_N = n-i,B_one = 1;
22       PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&B_N,&v2v,v1+i,&B_one,S+i*(n+1),&B_one));
23     }
24   }
25   PetscFunctionReturn(0);
26 }
27 
28 /* if v2 is not present, correction is done in-place */
29 #undef __FUNCT__
30 #define __FUNCT__ "PCBDDCReuseSolversChangeInterior"
31 PetscErrorCode PCBDDCReuseSolversChangeInterior(PCBDDCReuseSolvers ctx, Vec v, Vec v2, PetscBool sol)
32 {
33   PetscScalar    *array;
34   PetscScalar    *array2;
35   PetscErrorCode ierr;
36 
37   PetscFunctionBegin;
38   if (!ctx->benign_n) PetscFunctionReturn(0);
39   if (v2) {
40     PetscInt nl;
41 
42     ierr = VecGetArrayRead(v,(const PetscScalar**)&array);CHKERRQ(ierr);
43     ierr = VecGetLocalSize(v2,&nl);CHKERRQ(ierr);
44     ierr = VecGetArray(v2,&array2);CHKERRQ(ierr);
45     ierr = PetscMemcpy(array2,array,nl*sizeof(PetscScalar));CHKERRQ(ierr);
46   } else {
47     ierr = VecGetArray(v,&array);CHKERRQ(ierr);
48     array2 = array;
49   }
50   if (!sol) { /* change rhs */
51     PetscInt n;
52     for (n=0;n<ctx->benign_n;n++) {
53       PetscScalar    sum = 0.;
54       const PetscInt *cols;
55       PetscInt       nz,i;
56 
57       ierr = ISGetLocalSize(ctx->benign_zerodiag_subs[n],&nz);CHKERRQ(ierr);
58       ierr = ISGetIndices(ctx->benign_zerodiag_subs[n],&cols);CHKERRQ(ierr);
59       for (i=0;i<nz-1;i++) sum += array[cols[i]];
60       sum = -sum/nz;
61       for (i=0;i<nz-1;i++) array2[cols[i]] += sum;
62       ctx->benign_save_vals[n] = array2[cols[nz-1]];
63       array2[cols[nz-1]] = sum;
64       ierr = ISRestoreIndices(ctx->benign_zerodiag_subs[n],&cols);CHKERRQ(ierr);
65     }
66   } else {
67     PetscInt n;
68     for (n=0;n<ctx->benign_n;n++) {
69       PetscScalar    sum = 0.;
70       const PetscInt *cols;
71       PetscInt       nz,i;
72 
73       ierr = ISGetLocalSize(ctx->benign_zerodiag_subs[n],&nz);CHKERRQ(ierr);
74       ierr = ISGetIndices(ctx->benign_zerodiag_subs[n],&cols);CHKERRQ(ierr);
75       for (i=0;i<nz-1;i++) sum += array[cols[i]];
76       sum = -sum/nz;
77       for (i=0;i<nz-1;i++) array2[cols[i]] += sum;
78       array2[cols[nz-1]] = ctx->benign_save_vals[n];
79       ierr = ISRestoreIndices(ctx->benign_zerodiag_subs[n],&cols);CHKERRQ(ierr);
80     }
81   }
82   if (v2) {
83     ierr = VecRestoreArrayRead(v,(const PetscScalar**)&array);CHKERRQ(ierr);
84     ierr = VecRestoreArray(v2,&array);CHKERRQ(ierr);
85   } else {
86     ierr = VecRestoreArray(v,&array);CHKERRQ(ierr);
87   }
88   PetscFunctionReturn(0);
89 }
90 
91 #undef __FUNCT__
92 #define __FUNCT__ "PCBDDCReuseSolvers_Solve_Private"
93 static PetscErrorCode PCBDDCReuseSolvers_Solve_Private(PC pc, Vec rhs, Vec sol, PetscBool transpose, PetscBool full)
94 {
95   PCBDDCReuseSolvers ctx;
96   PetscBool          copy = PETSC_FALSE;
97   PetscErrorCode     ierr;
98 
99   PetscFunctionBegin;
100   ierr = PCShellGetContext(pc,(void **)&ctx);CHKERRQ(ierr);
101   if (full) {
102 #if defined(PETSC_HAVE_MUMPS)
103     ierr = MatMumpsSetIcntl(ctx->F,26,-1);CHKERRQ(ierr);
104 #endif
105     copy = ctx->has_vertices;
106   } else { /* interior solver */
107 #if defined(PETSC_HAVE_MUMPS)
108     ierr = MatMumpsSetIcntl(ctx->F,26,0);CHKERRQ(ierr);
109 #endif
110 #if defined(PETSC_HAVE_MKL_PARDISO)
111     ierr = MatMkl_PardisoSetCntl(ctx->F,70,1);CHKERRQ(ierr);
112 #endif
113     copy = PETSC_TRUE;
114   }
115   /* copy rhs into factored matrix workspace */
116   if (copy) {
117     PetscInt    n;
118     PetscScalar *array,*array_solver;
119 
120     ierr = VecGetLocalSize(rhs,&n);CHKERRQ(ierr);
121     ierr = VecGetArrayRead(rhs,(const PetscScalar**)&array);CHKERRQ(ierr);
122     ierr = VecGetArray(ctx->rhs,&array_solver);CHKERRQ(ierr);
123     ierr = PetscMemcpy(array_solver,array,n*sizeof(PetscScalar));CHKERRQ(ierr);
124     ierr = VecRestoreArray(ctx->rhs,&array_solver);CHKERRQ(ierr);
125     ierr = VecRestoreArrayRead(rhs,(const PetscScalar**)&array);CHKERRQ(ierr);
126 
127     ierr = PCBDDCReuseSolversChangeInterior(ctx,ctx->rhs,NULL,PETSC_FALSE);CHKERRQ(ierr);
128     if (transpose) {
129       ierr = MatSolveTranspose(ctx->F,ctx->rhs,ctx->sol);CHKERRQ(ierr);
130     } else {
131       ierr = MatSolve(ctx->F,ctx->rhs,ctx->sol);CHKERRQ(ierr);
132     }
133     ierr = PCBDDCReuseSolversChangeInterior(ctx,ctx->sol,NULL,PETSC_TRUE);CHKERRQ(ierr);
134 
135     /* get back data to caller worskpace */
136     ierr = VecGetArrayRead(ctx->sol,(const PetscScalar**)&array_solver);CHKERRQ(ierr);
137     ierr = VecGetArray(sol,&array);CHKERRQ(ierr);
138     ierr = PetscMemcpy(array,array_solver,n*sizeof(PetscScalar));CHKERRQ(ierr);
139     ierr = VecRestoreArray(sol,&array);CHKERRQ(ierr);
140     ierr = VecRestoreArrayRead(ctx->sol,(const PetscScalar**)&array_solver);CHKERRQ(ierr);
141   } else {
142     if (ctx->benign_n) {
143       ierr = PCBDDCReuseSolversChangeInterior(ctx,rhs,ctx->rhs,PETSC_FALSE);CHKERRQ(ierr);
144       if (transpose) {
145         ierr = MatSolveTranspose(ctx->F,ctx->rhs,sol);CHKERRQ(ierr);
146       } else {
147         ierr = MatSolve(ctx->F,ctx->rhs,sol);CHKERRQ(ierr);
148       }
149       ierr = PCBDDCReuseSolversChangeInterior(ctx,sol,NULL,PETSC_TRUE);CHKERRQ(ierr);
150     } else {
151       if (transpose) {
152         ierr = MatSolveTranspose(ctx->F,rhs,sol);CHKERRQ(ierr);
153       } else {
154         ierr = MatSolve(ctx->F,rhs,sol);CHKERRQ(ierr);
155       }
156     }
157   }
158 #if defined(PETSC_HAVE_MKL_PARDISO)
159   ierr = MatMkl_PardisoSetCntl(ctx->F,70,0);CHKERRQ(ierr);
160 #endif
161   PetscFunctionReturn(0);
162 }
163 
164 #undef __FUNCT__
165 #define __FUNCT__ "PCBDDCReuseSolvers_Correction"
166 static PetscErrorCode PCBDDCReuseSolvers_Correction(PC pc, Vec rhs, Vec sol)
167 {
168   PetscErrorCode   ierr;
169 
170   PetscFunctionBegin;
171   ierr = PCBDDCReuseSolvers_Solve_Private(pc,rhs,sol,PETSC_FALSE,PETSC_TRUE);CHKERRQ(ierr);
172   PetscFunctionReturn(0);
173 }
174 
175 #undef __FUNCT__
176 #define __FUNCT__ "PCBDDCReuseSolvers_CorrectionTranspose"
177 static PetscErrorCode PCBDDCReuseSolvers_CorrectionTranspose(PC pc, Vec rhs, Vec sol)
178 {
179   PetscErrorCode   ierr;
180 
181   PetscFunctionBegin;
182   ierr = PCBDDCReuseSolvers_Solve_Private(pc,rhs,sol,PETSC_TRUE,PETSC_TRUE);CHKERRQ(ierr);
183   PetscFunctionReturn(0);
184 }
185 
186 #undef __FUNCT__
187 #define __FUNCT__ "PCBDDCReuseSolvers_Interior"
188 static PetscErrorCode PCBDDCReuseSolvers_Interior(PC pc, Vec rhs, Vec sol)
189 {
190   PetscErrorCode   ierr;
191 
192   PetscFunctionBegin;
193   ierr = PCBDDCReuseSolvers_Solve_Private(pc,rhs,sol,PETSC_FALSE,PETSC_FALSE);CHKERRQ(ierr);
194   PetscFunctionReturn(0);
195 }
196 
197 #undef __FUNCT__
198 #define __FUNCT__ "PCBDDCReuseSolvers_InteriorTranspose"
199 static PetscErrorCode PCBDDCReuseSolvers_InteriorTranspose(PC pc, Vec rhs, Vec sol)
200 {
201   PetscErrorCode   ierr;
202 
203   PetscFunctionBegin;
204   ierr = PCBDDCReuseSolvers_Solve_Private(pc,rhs,sol,PETSC_TRUE,PETSC_FALSE);CHKERRQ(ierr);
205   PetscFunctionReturn(0);
206 }
207 
208 #undef __FUNCT__
209 #define __FUNCT__ "PCBDDCReuseSolversReset"
210 static PetscErrorCode PCBDDCReuseSolversReset(PCBDDCReuseSolvers reuse)
211 {
212   PetscInt       i;
213   PetscErrorCode ierr;
214 
215   PetscFunctionBegin;
216   ierr = MatDestroy(&reuse->F);CHKERRQ(ierr);
217   ierr = VecDestroy(&reuse->sol);CHKERRQ(ierr);
218   ierr = VecDestroy(&reuse->rhs);CHKERRQ(ierr);
219   ierr = PCDestroy(&reuse->interior_solver);CHKERRQ(ierr);
220   ierr = PCDestroy(&reuse->correction_solver);CHKERRQ(ierr);
221   ierr = ISDestroy(&reuse->is_R);CHKERRQ(ierr);
222   ierr = ISDestroy(&reuse->is_B);CHKERRQ(ierr);
223   ierr = VecScatterDestroy(&reuse->correction_scatter_B);CHKERRQ(ierr);
224   ierr = VecDestroy(&reuse->sol_B);CHKERRQ(ierr);
225   ierr = VecDestroy(&reuse->rhs_B);CHKERRQ(ierr);
226   for (i=0;i<reuse->benign_n;i++) {
227     ierr = ISDestroy(&reuse->benign_zerodiag_subs[i]);CHKERRQ(ierr);
228   }
229   ierr = PetscFree(reuse->benign_zerodiag_subs);CHKERRQ(ierr);
230   ierr = PetscFree(reuse->benign_save_vals);CHKERRQ(ierr);
231   PetscFunctionReturn(0);
232 }
233 
234 #undef __FUNCT__
235 #define __FUNCT__ "PCBDDCComputeExplicitSchur"
236 static PetscErrorCode PCBDDCComputeExplicitSchur(Mat M, PetscBool issym, MatReuse reuse, Mat *S)
237 {
238   Mat            B, C, D, Bd, Cd, AinvBd;
239   KSP            ksp;
240   PC             pc;
241   PetscBool      isLU, isILU, isCHOL, Bdense, Cdense;
242   PetscReal      fill = 2.0;
243   PetscInt       n_I;
244   PetscMPIInt    size;
245   PetscErrorCode ierr;
246 
247   PetscFunctionBegin;
248   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)M),&size);CHKERRQ(ierr);
249   if (size != 1) {
250     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not for parallel matrices");
251   }
252   if (reuse == MAT_REUSE_MATRIX) {
253     PetscBool Sdense;
254 
255     ierr = PetscObjectTypeCompare((PetscObject)*S, MATSEQDENSE, &Sdense);CHKERRQ(ierr);
256     if (!Sdense) {
257       SETERRQ(PetscObjectComm((PetscObject)M),PETSC_ERR_SUP,"S should be dense");
258     }
259   }
260   ierr = MatSchurComplementGetSubMatrices(M, NULL, NULL, &B, &C, &D);CHKERRQ(ierr);
261   ierr = MatSchurComplementGetKSP(M, &ksp);CHKERRQ(ierr);
262   ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr);
263   ierr = PetscObjectTypeCompare((PetscObject) pc, PCLU, &isLU);CHKERRQ(ierr);
264   ierr = PetscObjectTypeCompare((PetscObject) pc, PCILU, &isILU);CHKERRQ(ierr);
265   ierr = PetscObjectTypeCompare((PetscObject) pc, PCCHOLESKY, &isCHOL);CHKERRQ(ierr);
266   ierr = PetscObjectTypeCompare((PetscObject) B, MATSEQDENSE, &Bdense);CHKERRQ(ierr);
267   ierr = PetscObjectTypeCompare((PetscObject) C, MATSEQDENSE, &Cdense);CHKERRQ(ierr);
268   ierr = MatGetSize(B,&n_I,NULL);CHKERRQ(ierr);
269   if (n_I) {
270     if (!Bdense) {
271       ierr = MatConvert(B, MATSEQDENSE, MAT_INITIAL_MATRIX, &Bd);CHKERRQ(ierr);
272     } else {
273       Bd = B;
274     }
275 
276     if (isLU || isILU || isCHOL) {
277       Mat fact;
278       ierr = KSPSetUp(ksp);CHKERRQ(ierr);
279       ierr = PCFactorGetMatrix(pc, &fact);CHKERRQ(ierr);
280       ierr = MatDuplicate(Bd, MAT_DO_NOT_COPY_VALUES, &AinvBd);CHKERRQ(ierr);
281       ierr = MatMatSolve(fact, Bd, AinvBd);CHKERRQ(ierr);
282     } else {
283       PetscBool ex = PETSC_TRUE;
284 
285       if (ex) {
286         Mat Ainvd;
287 
288         ierr = PCComputeExplicitOperator(pc, &Ainvd);CHKERRQ(ierr);
289         ierr = MatMatMult(Ainvd, Bd, MAT_INITIAL_MATRIX, fill, &AinvBd);CHKERRQ(ierr);
290         ierr = MatDestroy(&Ainvd);CHKERRQ(ierr);
291       } else {
292         Vec         sol,rhs;
293         PetscScalar *arrayrhs,*arraysol;
294         PetscInt    i,nrhs,n;
295 
296         ierr = MatDuplicate(Bd, MAT_DO_NOT_COPY_VALUES, &AinvBd);CHKERRQ(ierr);
297         ierr = MatGetSize(Bd,&n,&nrhs);CHKERRQ(ierr);
298         ierr = MatDenseGetArray(Bd,&arrayrhs);CHKERRQ(ierr);
299         ierr = MatDenseGetArray(AinvBd,&arraysol);CHKERRQ(ierr);
300         ierr = KSPGetSolution(ksp,&sol);CHKERRQ(ierr);
301         ierr = KSPGetRhs(ksp,&rhs);CHKERRQ(ierr);
302         for (i=0;i<nrhs;i++) {
303           ierr = VecPlaceArray(rhs,arrayrhs+i*n);CHKERRQ(ierr);
304           ierr = VecPlaceArray(sol,arraysol+i*n);CHKERRQ(ierr);
305           ierr = KSPSolve(ksp,rhs,sol);CHKERRQ(ierr);
306           ierr = VecResetArray(rhs);CHKERRQ(ierr);
307           ierr = VecResetArray(sol);CHKERRQ(ierr);
308         }
309         ierr = MatDenseRestoreArray(Bd,&arrayrhs);CHKERRQ(ierr);
310         ierr = MatDenseRestoreArray(AinvBd,&arrayrhs);CHKERRQ(ierr);
311       }
312     }
313     if (!Bdense & !issym) {
314       ierr = MatDestroy(&Bd);CHKERRQ(ierr);
315     }
316 
317     if (!issym) {
318       if (!Cdense) {
319         ierr = MatConvert(C, MATSEQDENSE, MAT_INITIAL_MATRIX, &Cd);CHKERRQ(ierr);
320       } else {
321         Cd = C;
322       }
323       ierr = MatMatMult(Cd, AinvBd, reuse, fill, S);CHKERRQ(ierr);
324       if (!Cdense) {
325         ierr = MatDestroy(&Cd);CHKERRQ(ierr);
326       }
327     } else {
328       ierr = MatTransposeMatMult(Bd, AinvBd, reuse, fill, S);CHKERRQ(ierr);
329       if (!Bdense) {
330         ierr = MatDestroy(&Bd);CHKERRQ(ierr);
331       }
332     }
333     ierr = MatDestroy(&AinvBd);CHKERRQ(ierr);
334   }
335 
336   if (D) {
337     Mat       Dd;
338     PetscBool Ddense;
339 
340     ierr = PetscObjectTypeCompare((PetscObject)D,MATSEQDENSE,&Ddense);CHKERRQ(ierr);
341     if (!Ddense) {
342       ierr = MatConvert(D, MATSEQDENSE, MAT_INITIAL_MATRIX, &Dd);CHKERRQ(ierr);
343     } else {
344       Dd = D;
345     }
346     if (n_I) {
347       ierr = MatAYPX(*S,-1.0,Dd,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
348     } else {
349       if (reuse == MAT_INITIAL_MATRIX) {
350         ierr = MatDuplicate(Dd,MAT_COPY_VALUES,S);CHKERRQ(ierr);
351       } else {
352         ierr = MatCopy(Dd,*S,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
353       }
354     }
355     if (!Ddense) {
356       ierr = MatDestroy(&Dd);CHKERRQ(ierr);
357     }
358   } else {
359     ierr = MatScale(*S,-1.0);CHKERRQ(ierr);
360   }
361   PetscFunctionReturn(0);
362 }
363 
364 #undef __FUNCT__
365 #define __FUNCT__ "PCBDDCSubSchursSetUp"
366 PetscErrorCode PCBDDCSubSchursSetUp(PCBDDCSubSchurs sub_schurs, Mat Ain, Mat Sin, PetscBool exact_schur, PetscInt xadj[], PetscInt adjncy[], PetscInt nlayers, PetscBool faster_deluxe, PetscBool compute_Stilda, PetscBool reuse_solvers, PetscBool benign_trick, PetscInt benign_n, PetscInt benign_p0_lidx[], IS benign_zerodiag_subs[])
367 {
368   Mat                    F,A_II,A_IB,A_BI,A_BB,AE_II;
369   Mat                    S_all;
370   Mat                    global_schur_subsets,work_mat,*submats;
371   ISLocalToGlobalMapping l2gmap_subsets;
372   IS                     is_I,is_I_layer;
373   IS                     all_subsets,all_subsets_mult,all_subsets_n;
374   PetscInt               *nnz,*all_local_idx_N;
375   PetscInt               *auxnum1,*auxnum2;
376   PetscInt               i,subset_size,max_subset_size;
377   PetscInt               n_B,extra,local_size,global_size;
378   PetscBLASInt           B_N,B_ierr,B_lwork,*pivots;
379   PetscScalar            *Bwork;
380   PetscSubcomm           subcomm;
381   PetscMPIInt            color,rank;
382   MPI_Comm               comm_n;
383   PetscBool              use_getr = PETSC_FALSE;
384   PetscErrorCode         ierr;
385 
386   PetscFunctionBegin;
387   /* update info in sub_schurs */
388   ierr = MatDestroy(&sub_schurs->A);CHKERRQ(ierr);
389   ierr = MatDestroy(&sub_schurs->S);CHKERRQ(ierr);
390   sub_schurs->is_hermitian = PETSC_FALSE;
391   sub_schurs->is_posdef = PETSC_FALSE;
392   if (Ain) {
393     PetscBool isseqaij;
394     /* determine if we are dealing with hermitian positive definite problems */
395 #if !defined(PETSC_USE_COMPLEX)
396     if (Ain->symmetric_set) {
397       sub_schurs->is_hermitian = Ain->symmetric;
398     }
399 #else
400     if (Ain->hermitian_set) {
401       sub_schurs->is_hermitian = Ain->hermitian;
402     }
403 #endif
404     if (Ain->spd_set) {
405       sub_schurs->is_posdef = Ain->spd;
406     }
407 
408     /* check */
409     ierr = PetscObjectTypeCompare((PetscObject)Ain,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
410     if (compute_Stilda && (!Ain->hermitian_set || !Ain->spd_set)) { /* these are lazy checks, maybe I should throw an error if they are not set */
411       PetscInt lsize;
412 
413       ierr = MatGetSize(Ain,&lsize,NULL);CHKERRQ(ierr);
414       if (lsize) {
415         PetscScalar val;
416         PetscReal   norm;
417         Vec         vec1,vec2,vec3;
418 
419         ierr = MatCreateVecs(Ain,&vec1,&vec2);CHKERRQ(ierr);
420         ierr = VecDuplicate(vec1,&vec3);CHKERRQ(ierr);
421         ierr = VecSetRandom(vec1,NULL);CHKERRQ(ierr);
422         ierr = MatMult(Ain,vec1,vec2);CHKERRQ(ierr);
423 #if !defined(PETSC_USE_COMPLEX)
424         ierr = MatMultTranspose(Ain,vec1,vec3);CHKERRQ(ierr);
425 #else
426         ierr = MatMultHermitianTranspose(Ain,vec1,vec3);CHKERRQ(ierr);
427 #endif
428         ierr = VecAXPY(vec3,-1.0,vec2);CHKERRQ(ierr);
429         ierr = VecNorm(vec3,NORM_INFINITY,&norm);CHKERRQ(ierr);
430         if (!Ain->hermitian_set) {
431           if (norm > PetscSqrtReal(PETSC_SMALL)) {
432             sub_schurs->is_hermitian = PETSC_FALSE;
433           } else {
434             sub_schurs->is_hermitian = PETSC_TRUE;
435           }
436         }
437         if (!Ain->spd_set && !benign_trick) {
438           ierr = VecDot(vec1,vec2,&val);CHKERRQ(ierr);
439           if (PetscRealPart(val) > 0. && PetscAbsReal(PetscImaginaryPart(val)) < PETSC_SMALL) sub_schurs->is_posdef = PETSC_TRUE;
440         }
441         ierr = VecDestroy(&vec1);CHKERRQ(ierr);
442         ierr = VecDestroy(&vec2);CHKERRQ(ierr);
443         ierr = VecDestroy(&vec3);CHKERRQ(ierr);
444       } else {
445         sub_schurs->is_hermitian = PETSC_TRUE;
446         sub_schurs->is_posdef = PETSC_TRUE;
447       }
448     }
449     if (isseqaij) {
450       ierr = PetscObjectReference((PetscObject)Ain);CHKERRQ(ierr);
451       sub_schurs->A = Ain;
452     } else {
453       ierr = MatConvert(Ain,MATSEQAIJ,MAT_INITIAL_MATRIX,&sub_schurs->A);CHKERRQ(ierr);
454     }
455   }
456 
457   ierr = PetscObjectReference((PetscObject)Sin);CHKERRQ(ierr);
458   sub_schurs->S = Sin;
459   if (sub_schurs->schur_explicit) {
460     sub_schurs->schur_explicit = (PetscBool)(!!sub_schurs->A);
461   }
462 
463   /* preliminary checks */
464   if (!sub_schurs->schur_explicit && compute_Stilda) {
465     SETERRQ(PetscObjectComm((PetscObject)sub_schurs->l2gmap),PETSC_ERR_SUP,"Adaptive selection of constraints requires MUMPS and/or MKL_PARDISO");
466   }
467 
468   /* restrict work on active processes */
469   color = 0;
470   if (!sub_schurs->n_subs) color = 1; /* this can happen if we are in a multilevel case or if the subdomain is disconnected */
471   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&rank);CHKERRQ(ierr);
472   ierr = PetscSubcommCreate(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&subcomm);CHKERRQ(ierr);
473   ierr = PetscSubcommSetNumber(subcomm,2);CHKERRQ(ierr);
474   ierr = PetscSubcommSetTypeGeneral(subcomm,color,rank);CHKERRQ(ierr);
475   ierr = PetscCommDuplicate(PetscSubcommChild(subcomm),&comm_n,NULL);CHKERRQ(ierr);
476   ierr = PetscSubcommDestroy(&subcomm);CHKERRQ(ierr);
477   if (!sub_schurs->n_subs) {
478     ierr = PetscCommDestroy(&comm_n);CHKERRQ(ierr);
479     PetscFunctionReturn(0);
480   }
481   /* ierr = PetscCommDuplicate(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&comm_n,NULL);CHKERRQ(ierr); */
482 
483   /* get Schur complement matrices */
484   if (!sub_schurs->schur_explicit) {
485     Mat       tA_IB,tA_BI,tA_BB;
486     PetscBool isseqsbaij;
487     ierr = MatSchurComplementGetSubMatrices(sub_schurs->S,&A_II,NULL,&tA_IB,&tA_BI,&tA_BB);CHKERRQ(ierr);
488     ierr = PetscObjectTypeCompare((PetscObject)tA_BB,MATSEQSBAIJ,&isseqsbaij);CHKERRQ(ierr);
489     if (isseqsbaij) {
490       ierr = MatConvert(tA_BB,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_BB);CHKERRQ(ierr);
491       ierr = MatConvert(tA_IB,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_IB);CHKERRQ(ierr);
492       ierr = MatConvert(tA_BI,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_BI);CHKERRQ(ierr);
493     } else {
494       ierr = PetscObjectReference((PetscObject)tA_BB);CHKERRQ(ierr);
495       A_BB = tA_BB;
496       ierr = PetscObjectReference((PetscObject)tA_IB);CHKERRQ(ierr);
497       A_IB = tA_IB;
498       ierr = PetscObjectReference((PetscObject)tA_BI);CHKERRQ(ierr);
499       A_BI = tA_BI;
500     }
501   } else {
502     A_II = NULL;
503     A_IB = NULL;
504     A_BI = NULL;
505     A_BB = NULL;
506   }
507   S_all = NULL;
508 
509   /* determine interior problems */
510   ierr = ISGetLocalSize(sub_schurs->is_I,&i);CHKERRQ(ierr);
511   if (nlayers >= 0 && i) { /* Interior problems can be different from the original one */
512     PetscBT                touched;
513     const PetscInt*        idx_B;
514     PetscInt               n_I,n_B,n_local_dofs,n_prev_added,j,layer,*local_numbering;
515 
516     if (xadj == NULL || adjncy == NULL) {
517       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cannot request layering without adjacency");
518     }
519     /* get sizes */
520     ierr = ISGetLocalSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr);
521     ierr = ISGetLocalSize(sub_schurs->is_B,&n_B);CHKERRQ(ierr);
522 
523     ierr = PetscMalloc1(n_I+n_B,&local_numbering);CHKERRQ(ierr);
524     ierr = PetscBTCreate(n_I+n_B,&touched);CHKERRQ(ierr);
525     ierr = PetscBTMemzero(n_I+n_B,touched);CHKERRQ(ierr);
526 
527     /* all boundary dofs must be skipped when adding layers */
528     ierr = ISGetIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr);
529     for (j=0;j<n_B;j++) {
530       ierr = PetscBTSet(touched,idx_B[j]);CHKERRQ(ierr);
531     }
532     ierr = PetscMemcpy(local_numbering,idx_B,n_B*sizeof(PetscInt));CHKERRQ(ierr);
533     ierr = ISRestoreIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr);
534 
535     /* add prescribed number of layers of dofs */
536     n_local_dofs = n_B;
537     n_prev_added = n_B;
538     for (layer=0;layer<nlayers;layer++) {
539       PetscInt n_added;
540       if (n_local_dofs == n_I+n_B) break;
541       if (n_local_dofs > n_I+n_B) {
542         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);
543       }
544       ierr = PCBDDCAdjGetNextLayer_Private(local_numbering+n_local_dofs,n_prev_added,touched,xadj,adjncy,&n_added);CHKERRQ(ierr);
545       n_prev_added = n_added;
546       n_local_dofs += n_added;
547       if (!n_added) break;
548     }
549     ierr = PetscBTDestroy(&touched);CHKERRQ(ierr);
550 
551     /* IS for I layer dofs in original numbering */
552     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);
553     ierr = PetscFree(local_numbering);CHKERRQ(ierr);
554     ierr = ISSort(is_I_layer);CHKERRQ(ierr);
555     /* IS for I layer dofs in I numbering */
556     if (!sub_schurs->schur_explicit) {
557       ISLocalToGlobalMapping ItoNmap;
558       ierr = ISLocalToGlobalMappingCreateIS(sub_schurs->is_I,&ItoNmap);CHKERRQ(ierr);
559       ierr = ISGlobalToLocalMappingApplyIS(ItoNmap,IS_GTOLM_DROP,is_I_layer,&is_I);CHKERRQ(ierr);
560       ierr = ISLocalToGlobalMappingDestroy(&ItoNmap);CHKERRQ(ierr);
561 
562       /* II block */
563       ierr = MatGetSubMatrix(A_II,is_I,is_I,MAT_INITIAL_MATRIX,&AE_II);CHKERRQ(ierr);
564     }
565   } else {
566     PetscInt n_I;
567 
568     /* IS for I dofs in original numbering */
569     ierr = PetscObjectReference((PetscObject)sub_schurs->is_I);CHKERRQ(ierr);
570     is_I_layer = sub_schurs->is_I;
571 
572     /* IS for I dofs in I numbering (strided 1) */
573     if (!sub_schurs->schur_explicit) {
574       ierr = ISGetSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr);
575       ierr = ISCreateStride(PetscObjectComm((PetscObject)sub_schurs->is_I),n_I,0,1,&is_I);CHKERRQ(ierr);
576 
577       /* II block is the same */
578       ierr = PetscObjectReference((PetscObject)A_II);CHKERRQ(ierr);
579       AE_II = A_II;
580     }
581   }
582 
583   /* Get info on subset sizes and sum of all subsets sizes */
584   max_subset_size = 0;
585   local_size = 0;
586   for (i=0;i<sub_schurs->n_subs;i++) {
587     ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
588     max_subset_size = PetscMax(subset_size,max_subset_size);
589     local_size += subset_size;
590   }
591 
592   /* Work arrays for local indices */
593   extra = 0;
594   ierr = ISGetLocalSize(sub_schurs->is_B,&n_B);CHKERRQ(ierr);
595   if (sub_schurs->schur_explicit && is_I_layer) {
596     ierr = ISGetLocalSize(is_I_layer,&extra);CHKERRQ(ierr);
597   }
598   ierr = PetscMalloc1(n_B+extra,&all_local_idx_N);CHKERRQ(ierr);
599   if (extra) {
600     const PetscInt *idxs;
601     ierr = ISGetIndices(is_I_layer,&idxs);CHKERRQ(ierr);
602     ierr = PetscMemcpy(all_local_idx_N,idxs,extra*sizeof(PetscInt));CHKERRQ(ierr);
603     ierr = ISRestoreIndices(is_I_layer,&idxs);CHKERRQ(ierr);
604   }
605   ierr = PetscMalloc1(local_size,&nnz);CHKERRQ(ierr);
606   ierr = PetscMalloc1(sub_schurs->n_subs,&auxnum1);CHKERRQ(ierr);
607   ierr = PetscMalloc1(sub_schurs->n_subs,&auxnum2);CHKERRQ(ierr);
608 
609   /* Get local indices in local numbering */
610   local_size = 0;
611   for (i=0;i<sub_schurs->n_subs;i++) {
612     PetscInt j;
613     const    PetscInt *idxs;
614 
615     ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
616     ierr = ISGetIndices(sub_schurs->is_subs[i],&idxs);CHKERRQ(ierr);
617     /* start (smallest in global ordering) and multiplicity */
618     auxnum1[i] = idxs[0];
619     auxnum2[i] = subset_size;
620     /* subset indices in local numbering */
621     ierr = PetscMemcpy(all_local_idx_N+local_size+extra,idxs,subset_size*sizeof(PetscInt));CHKERRQ(ierr);
622     ierr = ISRestoreIndices(sub_schurs->is_subs[i],&idxs);CHKERRQ(ierr);
623     for (j=0;j<subset_size;j++) nnz[local_size+j] = subset_size;
624     local_size += subset_size;
625   }
626 
627   /* allocate extra workspace needed only for GETRI */
628   Bwork = NULL;
629   pivots = NULL;
630   if (local_size && !benign_trick && (!sub_schurs->is_hermitian || !sub_schurs->is_posdef)) {
631     PetscScalar lwork;
632 
633     use_getr = PETSC_TRUE;
634     B_lwork = -1;
635     ierr = PetscBLASIntCast(local_size,&B_N);CHKERRQ(ierr);
636     ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
637     PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,Bwork,&B_N,pivots,&lwork,&B_lwork,&B_ierr));
638     ierr = PetscFPTrapPop();CHKERRQ(ierr);
639     if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to GETRI Lapack routine %d",(int)B_ierr);
640     ierr = PetscBLASIntCast((PetscInt)PetscRealPart(lwork),&B_lwork);CHKERRQ(ierr);
641     ierr = PetscMalloc2(B_lwork,&Bwork,B_N,&pivots);CHKERRQ(ierr);
642   }
643 
644   /* prepare parallel matrices for summing up properly schurs on subsets */
645   ierr = ISCreateGeneral(comm_n,sub_schurs->n_subs,auxnum1,PETSC_OWN_POINTER,&all_subsets_n);CHKERRQ(ierr);
646   ierr = ISLocalToGlobalMappingApplyIS(sub_schurs->l2gmap,all_subsets_n,&all_subsets);CHKERRQ(ierr);
647   ierr = ISDestroy(&all_subsets_n);CHKERRQ(ierr);
648   ierr = ISCreateGeneral(comm_n,sub_schurs->n_subs,auxnum2,PETSC_OWN_POINTER,&all_subsets_mult);CHKERRQ(ierr);
649   ierr = PCBDDCSubsetNumbering(all_subsets,all_subsets_mult,&global_size,&all_subsets_n);CHKERRQ(ierr);
650   ierr = ISDestroy(&all_subsets);CHKERRQ(ierr);
651   ierr = ISDestroy(&all_subsets_mult);CHKERRQ(ierr);
652   ierr = ISGetLocalSize(all_subsets_n,&i);CHKERRQ(ierr);
653   if (i != local_size) {
654     SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Invalid size of new subset! %d != %d",i,local_size);
655   }
656   ierr = ISLocalToGlobalMappingCreateIS(all_subsets_n,&l2gmap_subsets);CHKERRQ(ierr);
657   ierr = MatCreateIS(comm_n,1,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size,l2gmap_subsets,NULL,&work_mat);CHKERRQ(ierr);
658   ierr = ISLocalToGlobalMappingDestroy(&l2gmap_subsets);CHKERRQ(ierr);
659   ierr = MatCreate(PetscObjectComm((PetscObject)work_mat),&global_schur_subsets);CHKERRQ(ierr);
660   ierr = MatSetSizes(global_schur_subsets,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size);CHKERRQ(ierr);
661   ierr = MatSetType(global_schur_subsets,MATMPIAIJ);CHKERRQ(ierr);
662 
663   /* subset indices in local boundary numbering */
664   if (!sub_schurs->is_Ej_all) {
665     PetscInt *all_local_idx_B;
666 
667     ierr = PetscMalloc1(local_size,&all_local_idx_B);CHKERRQ(ierr);
668     ierr = ISGlobalToLocalMappingApply(sub_schurs->BtoNmap,IS_GTOLM_DROP,local_size,all_local_idx_N+extra,&subset_size,all_local_idx_B);CHKERRQ(ierr);
669     if (subset_size != local_size) {
670       SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in sub_schurs serial (BtoNmap)! %d != %d\n",subset_size,local_size);
671     }
672     ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size,all_local_idx_B,PETSC_OWN_POINTER,&sub_schurs->is_Ej_all);CHKERRQ(ierr);
673   }
674 
675   /* Local matrix of all local Schur on subsets (transposed) */
676   if (!sub_schurs->S_Ej_all) {
677     ierr = MatCreate(PETSC_COMM_SELF,&sub_schurs->S_Ej_all);CHKERRQ(ierr);
678     ierr = MatSetSizes(sub_schurs->S_Ej_all,PETSC_DECIDE,PETSC_DECIDE,local_size,local_size);CHKERRQ(ierr);
679     ierr = MatSetType(sub_schurs->S_Ej_all,MATAIJ);CHKERRQ(ierr);
680     ierr = MatSeqAIJSetPreallocation(sub_schurs->S_Ej_all,0,nnz);CHKERRQ(ierr);
681   }
682 
683   /* Compute Schur complements explicitly */
684   F = NULL;
685   if (!sub_schurs->schur_explicit) { /* this code branch is used when MatFactor with Schur complemnent support is not present; it is not very efficient, unless the economic version of the scaling is required */
686     Mat         S_Ej_expl;
687     PetscScalar *work;
688     PetscInt    j,*dummy_idx;
689     PetscBool   Sdense;
690 
691     ierr = PetscMalloc2(max_subset_size,&dummy_idx,max_subset_size*max_subset_size,&work);CHKERRQ(ierr);
692     local_size = 0;
693     for (i=0;i<sub_schurs->n_subs;i++) {
694       IS  is_subset_B;
695       Mat AE_EE,AE_IE,AE_EI,S_Ej;
696 
697       /* subsets in original and boundary numbering */
698       ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,sub_schurs->is_subs[i],&is_subset_B);CHKERRQ(ierr);
699       /* EE block */
700       ierr = MatGetSubMatrix(A_BB,is_subset_B,is_subset_B,MAT_INITIAL_MATRIX,&AE_EE);CHKERRQ(ierr);
701       /* IE block */
702       ierr = MatGetSubMatrix(A_IB,is_I,is_subset_B,MAT_INITIAL_MATRIX,&AE_IE);CHKERRQ(ierr);
703       /* EI block */
704       if (sub_schurs->is_hermitian) {
705         ierr = MatCreateTranspose(AE_IE,&AE_EI);CHKERRQ(ierr);
706       } else {
707         ierr = MatGetSubMatrix(A_BI,is_subset_B,is_I,MAT_INITIAL_MATRIX,&AE_EI);CHKERRQ(ierr);
708       }
709       ierr = ISDestroy(&is_subset_B);CHKERRQ(ierr);
710       ierr = MatCreateSchurComplement(AE_II,AE_II,AE_IE,AE_EI,AE_EE,&S_Ej);CHKERRQ(ierr);
711       ierr = MatDestroy(&AE_EE);CHKERRQ(ierr);
712       ierr = MatDestroy(&AE_IE);CHKERRQ(ierr);
713       ierr = MatDestroy(&AE_EI);CHKERRQ(ierr);
714       if (AE_II == A_II) { /* we can reuse the same ksp */
715         KSP ksp;
716         ierr = MatSchurComplementGetKSP(sub_schurs->S,&ksp);CHKERRQ(ierr);
717         ierr = MatSchurComplementSetKSP(S_Ej,ksp);CHKERRQ(ierr);
718       } else { /* build new ksp object which inherits ksp and pc types from the original one */
719         KSP       origksp,schurksp;
720         PC        origpc,schurpc;
721         KSPType   ksp_type;
722         PetscInt  n_internal;
723         PetscBool ispcnone;
724 
725         ierr = MatSchurComplementGetKSP(sub_schurs->S,&origksp);CHKERRQ(ierr);
726         ierr = MatSchurComplementGetKSP(S_Ej,&schurksp);CHKERRQ(ierr);
727         ierr = KSPGetType(origksp,&ksp_type);CHKERRQ(ierr);
728         ierr = KSPSetType(schurksp,ksp_type);CHKERRQ(ierr);
729         ierr = KSPGetPC(schurksp,&schurpc);CHKERRQ(ierr);
730         ierr = KSPGetPC(origksp,&origpc);CHKERRQ(ierr);
731         ierr = PetscObjectTypeCompare((PetscObject)origpc,PCNONE,&ispcnone);CHKERRQ(ierr);
732         if (!ispcnone) {
733           PCType pc_type;
734           ierr = PCGetType(origpc,&pc_type);CHKERRQ(ierr);
735           ierr = PCSetType(schurpc,pc_type);CHKERRQ(ierr);
736         } else {
737           ierr = PCSetType(schurpc,PCLU);CHKERRQ(ierr);
738         }
739         ierr = ISGetSize(is_I,&n_internal);CHKERRQ(ierr);
740         if (n_internal) { /* UMFPACK gives error with 0 sized problems */
741           MatSolverPackage solver=NULL;
742           ierr = PCFactorGetMatSolverPackage(origpc,(const MatSolverPackage*)&solver);CHKERRQ(ierr);
743           if (solver) {
744             ierr = PCFactorSetMatSolverPackage(schurpc,solver);CHKERRQ(ierr);
745           }
746         }
747         ierr = KSPSetUp(schurksp);CHKERRQ(ierr);
748       }
749       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
750       ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&S_Ej_expl);CHKERRQ(ierr);
751       ierr = PCBDDCComputeExplicitSchur(S_Ej,sub_schurs->is_hermitian,MAT_REUSE_MATRIX,&S_Ej_expl);CHKERRQ(ierr);
752       ierr = PetscObjectTypeCompare((PetscObject)S_Ej_expl,MATSEQDENSE,&Sdense);CHKERRQ(ierr);
753       if (Sdense) {
754         for (j=0;j<subset_size;j++) {
755           dummy_idx[j]=local_size+j;
756         }
757         ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
758       } else {
759         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not yet implemented for sparse matrices");
760       }
761       ierr = MatDestroy(&S_Ej);CHKERRQ(ierr);
762       ierr = MatDestroy(&S_Ej_expl);CHKERRQ(ierr);
763       local_size += subset_size;
764     }
765     ierr = PetscFree2(dummy_idx,work);CHKERRQ(ierr);
766     /* free */
767     ierr = ISDestroy(&is_I);CHKERRQ(ierr);
768     ierr = MatDestroy(&AE_II);CHKERRQ(ierr);
769     ierr = PetscFree(all_local_idx_N);CHKERRQ(ierr);
770   } else {
771     Mat         A;
772     Vec         *benign_AIIm1_ones = NULL;
773     IS          is_A_all,*is_p_r = NULL;
774     PetscScalar *work,*S_data,*schur_factor,*cs_AIB = NULL;
775     PetscInt    n,n_I,*dummy_idx,size_schur,size_active_schur,cum,cum2;
776     PetscBool   economic,solver_S,S_lower_triangular = PETSC_FALSE,factor_workaround;
777     char        solver[256];
778 
779     /* get sizes */
780     n_I = 0;
781     if (is_I_layer) {
782       ierr = ISGetLocalSize(is_I_layer,&n_I);CHKERRQ(ierr);
783     }
784     economic = PETSC_FALSE;
785     ierr = ISGetLocalSize(sub_schurs->is_I,&cum);CHKERRQ(ierr);
786     if (cum != n_I) economic = PETSC_TRUE;
787     ierr = MatGetLocalSize(sub_schurs->A,&n,NULL);CHKERRQ(ierr);
788 
789     /* size active schurs does not count any dirichlet or vertex dof on the interface */
790     size_active_schur = local_size;
791     cum = n_I+size_active_schur;
792     if (sub_schurs->is_dir) {
793       const PetscInt* idxs;
794       PetscInt        n_dir;
795 
796       ierr = ISGetLocalSize(sub_schurs->is_dir,&n_dir);CHKERRQ(ierr);
797       ierr = ISGetIndices(sub_schurs->is_dir,&idxs);CHKERRQ(ierr);
798       ierr = PetscMemcpy(all_local_idx_N+cum,idxs,n_dir*sizeof(PetscInt));CHKERRQ(ierr);
799       ierr = ISRestoreIndices(sub_schurs->is_dir,&idxs);CHKERRQ(ierr);
800       cum += n_dir;
801     }
802     factor_workaround = PETSC_FALSE;
803     /* include the primal vertices in the Schur complement */
804     if (exact_schur && sub_schurs->is_vertices && compute_Stilda) {
805       PetscInt n_v;
806 
807       ierr = ISGetLocalSize(sub_schurs->is_vertices,&n_v);CHKERRQ(ierr);
808       if (n_v) {
809         const PetscInt* idxs;
810 
811         ierr = ISGetIndices(sub_schurs->is_vertices,&idxs);CHKERRQ(ierr);
812         ierr = PetscMemcpy(all_local_idx_N+cum,idxs,n_v*sizeof(PetscInt));CHKERRQ(ierr);
813         ierr = ISRestoreIndices(sub_schurs->is_vertices,&idxs);CHKERRQ(ierr);
814         cum += n_v;
815         factor_workaround = PETSC_TRUE;
816       }
817     }
818     size_schur = cum - n_I;
819     ierr = ISCreateGeneral(PETSC_COMM_SELF,cum,all_local_idx_N,PETSC_OWN_POINTER,&is_A_all);CHKERRQ(ierr);
820     /* get working mat (TODO: factorize without actually permuting)  */
821     if (cum == n) {
822       ierr = ISSetPermutation(is_A_all);CHKERRQ(ierr);
823       ierr = MatPermute(sub_schurs->A,is_A_all,is_A_all,&A);CHKERRQ(ierr);
824     } else {
825       ierr = MatGetSubMatrix(sub_schurs->A,is_A_all,is_A_all,MAT_INITIAL_MATRIX,&A);CHKERRQ(ierr);
826     }
827     ierr = MatSetOptionsPrefix(A,"sub_schurs_");CHKERRQ(ierr);
828 
829     /* if we actually change the basis for the pressures, LDL^T factors will use a lot of memory
830        more than n^2 instead of n^1.5 or something. This is a workaround */
831     if (benign_n) {
832       Vec                    v;
833       ISLocalToGlobalMapping N_to_reor;
834       IS                     is_p0,is_p0_p;
835 
836       ierr = ISLocalToGlobalMappingCreateIS(is_A_all,&N_to_reor);CHKERRQ(ierr);
837       ierr = ISCreateGeneral(PETSC_COMM_SELF,benign_n,benign_p0_lidx,PETSC_COPY_VALUES,&is_p0);CHKERRQ(ierr);
838       ierr = ISGlobalToLocalMappingApplyIS(N_to_reor,IS_GTOLM_DROP,is_p0,&is_p0_p);CHKERRQ(ierr);
839       ierr = ISDestroy(&is_p0);CHKERRQ(ierr);
840       ierr = MatCreateVecs(A,&v,NULL);CHKERRQ(ierr);
841       ierr = VecDuplicateVecs(v,benign_n,&benign_AIIm1_ones);CHKERRQ(ierr);
842       ierr = PetscMalloc1(size_schur*benign_n,&cs_AIB);CHKERRQ(ierr);
843       ierr = PetscMalloc1(benign_n,&is_p_r);CHKERRQ(ierr);
844       /* compute colsum of A_IB restricted to pressures */
845       for (i=0;i<benign_n;i++) {
846         PetscScalar    *array;
847         const PetscInt *idxs;
848         PetscInt       j,nz;
849 
850         ierr = ISGlobalToLocalMappingApplyIS(N_to_reor,IS_GTOLM_DROP,benign_zerodiag_subs[i],&is_p_r[i]);CHKERRQ(ierr);
851         ierr = VecGetArray(benign_AIIm1_ones[i],&array);CHKERRQ(ierr);
852         ierr = ISGetLocalSize(is_p_r[i],&nz);CHKERRQ(ierr);
853         ierr = ISGetIndices(is_p_r[i],&idxs);CHKERRQ(ierr);
854         for (j=0;j<nz;j++) array[idxs[j]] = 1.;
855         ierr = ISRestoreIndices(is_p_r[i],&idxs);CHKERRQ(ierr);
856         ierr = VecRestoreArray(benign_AIIm1_ones[i],&array);CHKERRQ(ierr);
857         ierr = MatMult(A,benign_AIIm1_ones[i],v);CHKERRQ(ierr);
858         ierr = VecGetArray(v,&array);CHKERRQ(ierr);
859         for (j=0;j<size_schur;j++) cs_AIB[i*size_schur + j] = array[j+n_I];
860         ierr = VecRestoreArray(v,&array);CHKERRQ(ierr);
861       }
862       ierr = VecDestroy(&v);CHKERRQ(ierr);
863       ierr = MatSetOption(A,MAT_KEEP_NONZERO_PATTERN,PETSC_FALSE);CHKERRQ(ierr);
864       ierr = MatZeroRowsColumnsIS(A,is_p0_p,1.0,NULL,NULL);CHKERRQ(ierr);
865       ierr = ISDestroy(&is_p0_p);CHKERRQ(ierr);
866       ierr = ISLocalToGlobalMappingDestroy(&N_to_reor);CHKERRQ(ierr);
867     }
868     ierr = MatSetOption(A,MAT_SYMMETRIC,sub_schurs->is_hermitian);CHKERRQ(ierr);
869     ierr = MatSetOption(A,MAT_HERMITIAN,sub_schurs->is_hermitian);CHKERRQ(ierr);
870     ierr = MatSetOption(A,MAT_SPD,sub_schurs->is_posdef);CHKERRQ(ierr);
871 #if defined(PETSC_HAVE_MUMPS)
872     ierr = PetscStrcpy(solver,MATSOLVERMUMPS);CHKERRQ(ierr);
873 #else
874     ierr = PetscStrcpy(solver,MATSOLVERMKL_PARDISO);CHKERRQ(ierr);
875 #endif
876     ierr = PetscOptionsGetString("sub_schurs_","-mat_solver_package",solver,256,NULL);CHKERRQ(ierr);
877 
878     /* when using the benign subspace trick, the local Schur complements are SPD */
879     if (benign_trick) sub_schurs->is_posdef = PETSC_TRUE;
880 
881     if (n_I) { /* TODO add ordering from options */
882       IS is_schur;
883 
884       if (sub_schurs->is_hermitian && sub_schurs->is_posdef) {
885         ierr = MatGetFactor(A,solver,MAT_FACTOR_CHOLESKY,&F);CHKERRQ(ierr);
886       } else {
887         ierr = MatGetFactor(A,solver,MAT_FACTOR_LU,&F);CHKERRQ(ierr);
888       }
889       /* subsets ordered last */
890       ierr = ISCreateStride(PETSC_COMM_SELF,size_schur,n_I,1,&is_schur);CHKERRQ(ierr);
891       ierr = MatFactorSetSchurIS(F,is_schur);CHKERRQ(ierr);
892       ierr = ISDestroy(&is_schur);CHKERRQ(ierr);
893 
894       /* factorization step */
895       if (sub_schurs->is_hermitian && sub_schurs->is_posdef) {
896         ierr = MatCholeskyFactorSymbolic(F,A,NULL,NULL);CHKERRQ(ierr);
897 #if defined(PETSC_HAVE_MUMPS) /* be sure that icntl 19 is not set by command line */
898         ierr = MatMumpsSetIcntl(F,19,2);CHKERRQ(ierr);
899 #endif
900         if (sub_schurs->is_posdef) {
901           ierr = MatFactorSetSchurComplementSolverType(F,1);CHKERRQ(ierr);
902         } else {
903           ierr = MatFactorSetSchurComplementSolverType(F,2);CHKERRQ(ierr);
904         }
905         ierr = MatCholeskyFactorNumeric(F,A,NULL);CHKERRQ(ierr);
906         S_lower_triangular = PETSC_TRUE;
907       } else {
908         ierr = MatLUFactorSymbolic(F,A,NULL,NULL,NULL);CHKERRQ(ierr);
909 #if defined(PETSC_HAVE_MUMPS) /* be sure that icntl 19 is not set by command line */
910         ierr = MatMumpsSetIcntl(F,19,3);CHKERRQ(ierr);
911 #endif
912         ierr = MatLUFactorNumeric(F,A,NULL);CHKERRQ(ierr);
913       }
914 
915       /* get explicit Schur Complement computed during numeric factorization */
916       ierr = MatFactorGetSchurComplement(F,&S_all);CHKERRQ(ierr);
917       /* we can reuse the solvers if we are not using the economic version */
918       reuse_solvers = (PetscBool)(reuse_solvers && !economic);
919       factor_workaround = (PetscBool)(reuse_solvers && factor_workaround);
920       solver_S = PETSC_TRUE;
921 
922       /* update the Schur complement with the change of basis */
923       if (benign_n) {
924         PetscScalar *S_data;
925         Vec         v;
926 
927         ierr = MatDenseGetArray(S_all,&S_data);CHKERRQ(ierr);
928         ierr = VecDuplicate(benign_AIIm1_ones[0],&v);CHKERRQ(ierr);
929 #if defined(PETSC_HAVE_MUMPS)
930         ierr = MatMumpsSetIcntl(F,26,0);CHKERRQ(ierr);
931 #endif
932 #if defined(PETSC_HAVE_MKL_PARDISO)
933         ierr = MatMkl_PardisoSetCntl(F,70,1);CHKERRQ(ierr);
934 #endif
935         for (i=0;i<benign_n;i++) {
936           PetscScalar    *array,sum = 0.;
937           const PetscInt *idxs;
938           PetscInt       j,nz;
939 
940           ierr = VecCopy(benign_AIIm1_ones[i],v);CHKERRQ(ierr);
941           ierr = MatSolve(F,v,benign_AIIm1_ones[i]);CHKERRQ(ierr);
942           ierr = VecGetArray(benign_AIIm1_ones[i],&array);CHKERRQ(ierr);
943           ierr = ISGetLocalSize(is_p_r[i],&nz);CHKERRQ(ierr);
944           ierr = ISGetIndices(is_p_r[i],&idxs);CHKERRQ(ierr);
945           for (j=0;j<nz;j++) sum += array[idxs[j]];
946           sum -= 1.; /* p0 dof (eliminated) is excluded from the sum */
947           ierr = ISRestoreIndices(is_p_r[i],&idxs);CHKERRQ(ierr);
948           ierr = VecRestoreArray(benign_AIIm1_ones[i],&array);CHKERRQ(ierr);
949           ierr = MatMult(A,benign_AIIm1_ones[i],v);CHKERRQ(ierr);
950           ierr = VecGetArray(v,&array);CHKERRQ(ierr);
951           /* perform sparse rank-1 updates on symmetric Schur */
952           sum  = -sum/(1.*nz*nz);
953           ierr = SparseRankOneUpdate(S_data,size_schur,sum,cs_AIB+i*size_schur,cs_AIB+i*size_schur);CHKERRQ(ierr);
954           sum  = 1./nz;
955           ierr = SparseRankOneUpdate(S_data,size_schur,sum,array+n_I,cs_AIB+i*size_schur);CHKERRQ(ierr);
956           ierr = SparseRankOneUpdate(S_data,size_schur,sum,cs_AIB+i*size_schur,array+n_I);CHKERRQ(ierr);
957           ierr = VecRestoreArray(v,&array);CHKERRQ(ierr);
958         }
959         ierr = VecDestroy(&v);CHKERRQ(ierr);
960         ierr = MatDenseRestoreArray(S_all,&S_data);CHKERRQ(ierr);
961       }
962       ierr = PetscFree(cs_AIB);CHKERRQ(ierr);
963       ierr = VecDestroyVecs(benign_n,&benign_AIIm1_ones);CHKERRQ(ierr);
964       if (!reuse_solvers) {
965         for (i=0;i<benign_n;i++) {
966           ierr = ISDestroy(&is_p_r[i]);CHKERRQ(ierr);
967         }
968         ierr = PetscFree(is_p_r);CHKERRQ(ierr);
969       }
970     } else { /* we can't use MatFactor when size_schur == size_of_the_problem */
971       ierr = MatConvert(A,MATSEQDENSE,MAT_INITIAL_MATRIX,&S_all);CHKERRQ(ierr);
972       reuse_solvers = PETSC_FALSE; /* TODO: why we can't reuse the solvers here? */
973       solver_S = PETSC_FALSE;
974     }
975 
976     if (reuse_solvers) {
977       Mat                A_II;
978       Vec                vec1_B;
979       PCBDDCReuseSolvers msolv_ctx;
980 
981       if (sub_schurs->reuse_solver) {
982         ierr = PCBDDCReuseSolversReset(sub_schurs->reuse_solver);CHKERRQ(ierr);
983       } else {
984         ierr = PetscNew(&sub_schurs->reuse_solver);CHKERRQ(ierr);
985       }
986       msolv_ctx = sub_schurs->reuse_solver;
987       ierr = MatSchurComplementGetSubMatrices(sub_schurs->S,&A_II,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
988       ierr = PetscObjectReference((PetscObject)F);CHKERRQ(ierr);
989       msolv_ctx->F = F;
990       ierr = MatCreateVecs(F,&msolv_ctx->sol,NULL);CHKERRQ(ierr);
991       /* currently PETSc has no support for MatSolve(F,x,x), so cheat and let rhs and sol share the same memory */
992       {
993         PetscScalar *array;
994         PetscInt    n;
995 
996         ierr = VecGetLocalSize(msolv_ctx->sol,&n);CHKERRQ(ierr);
997         ierr = VecGetArray(msolv_ctx->sol,&array);CHKERRQ(ierr);
998         ierr = VecCreateSeqWithArray(PetscObjectComm((PetscObject)msolv_ctx->sol),1,n,array,&msolv_ctx->rhs);CHKERRQ(ierr);
999         ierr = VecRestoreArray(msolv_ctx->sol,&array);CHKERRQ(ierr);
1000       }
1001       msolv_ctx->has_vertices = factor_workaround;
1002 
1003       /* interior solver */
1004       ierr = PCCreate(PetscObjectComm((PetscObject)A_II),&msolv_ctx->interior_solver);CHKERRQ(ierr);
1005       ierr = PCSetOperators(msolv_ctx->interior_solver,A_II,A_II);CHKERRQ(ierr);
1006       ierr = PCSetType(msolv_ctx->interior_solver,PCSHELL);CHKERRQ(ierr);
1007       ierr = PCShellSetContext(msolv_ctx->interior_solver,msolv_ctx);CHKERRQ(ierr);
1008       ierr = PCShellSetApply(msolv_ctx->interior_solver,PCBDDCReuseSolvers_Interior);CHKERRQ(ierr);
1009       ierr = PCShellSetApplyTranspose(msolv_ctx->interior_solver,PCBDDCReuseSolvers_InteriorTranspose);CHKERRQ(ierr);
1010 
1011       /* correction solver */
1012       ierr = PCCreate(PetscObjectComm((PetscObject)A),&msolv_ctx->correction_solver);CHKERRQ(ierr);
1013       ierr = PCSetOperators(msolv_ctx->correction_solver,A,A);CHKERRQ(ierr);
1014       ierr = PCSetType(msolv_ctx->correction_solver,PCSHELL);CHKERRQ(ierr);
1015       ierr = PCShellSetContext(msolv_ctx->correction_solver,msolv_ctx);CHKERRQ(ierr);
1016       ierr = PCShellSetApply(msolv_ctx->correction_solver,PCBDDCReuseSolvers_Correction);CHKERRQ(ierr);
1017       ierr = PCShellSetApplyTranspose(msolv_ctx->correction_solver,PCBDDCReuseSolvers_CorrectionTranspose);CHKERRQ(ierr);
1018 
1019       /* scatter and vecs for Schur complement solver */
1020       ierr = MatCreateVecs(S_all,&msolv_ctx->sol_B,&msolv_ctx->rhs_B);CHKERRQ(ierr);
1021       ierr = MatCreateVecs(sub_schurs->S,&vec1_B,NULL);CHKERRQ(ierr);
1022       if (!factor_workaround) {
1023         ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,is_A_all,&msolv_ctx->is_B);CHKERRQ(ierr);
1024         ierr = VecScatterCreate(vec1_B,msolv_ctx->is_B,msolv_ctx->sol_B,NULL,&msolv_ctx->correction_scatter_B);CHKERRQ(ierr);
1025         ierr = PetscObjectReference((PetscObject)is_A_all);CHKERRQ(ierr);
1026         msolv_ctx->is_R = is_A_all;
1027       } else {
1028         IS              is_B_all;
1029         const PetscInt* idxs;
1030         PetscInt        dual,n_v,n;
1031 
1032         ierr = ISGetLocalSize(sub_schurs->is_vertices,&n_v);CHKERRQ(ierr);
1033         dual = size_schur - n_v;
1034         ierr = ISGetLocalSize(is_A_all,&n);CHKERRQ(ierr);
1035         ierr = ISGetIndices(is_A_all,&idxs);CHKERRQ(ierr);
1036         ierr = ISCreateGeneral(PetscObjectComm((PetscObject)is_A_all),dual,idxs+n_I,PETSC_COPY_VALUES,&is_B_all);CHKERRQ(ierr);
1037         ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,is_B_all,&msolv_ctx->is_B);CHKERRQ(ierr);
1038         ierr = ISDestroy(&is_B_all);CHKERRQ(ierr);
1039         ierr = ISCreateStride(PetscObjectComm((PetscObject)is_A_all),dual,0,1,&is_B_all);CHKERRQ(ierr);
1040         ierr = VecScatterCreate(vec1_B,msolv_ctx->is_B,msolv_ctx->sol_B,is_B_all,&msolv_ctx->correction_scatter_B);CHKERRQ(ierr);
1041         ierr = ISDestroy(&is_B_all);CHKERRQ(ierr);
1042         ierr = ISCreateGeneral(PetscObjectComm((PetscObject)is_A_all),n-n_v,idxs,PETSC_COPY_VALUES,&msolv_ctx->is_R);CHKERRQ(ierr);
1043         ierr = ISRestoreIndices(is_A_all,&idxs);CHKERRQ(ierr);
1044       }
1045       ierr = VecDestroy(&vec1_B);CHKERRQ(ierr);
1046 
1047       /* communicate benign info to solver context */
1048       if (benign_n) {
1049         msolv_ctx->benign_n = benign_n;
1050         msolv_ctx->benign_zerodiag_subs = is_p_r;
1051         ierr = PetscMalloc1(benign_n,&msolv_ctx->benign_save_vals);CHKERRQ(ierr);
1052       }
1053     }
1054     ierr = MatDestroy(&A);CHKERRQ(ierr);
1055     ierr = ISDestroy(&is_A_all);CHKERRQ(ierr);
1056 
1057     /* Work arrays */
1058     ierr = PetscMalloc2(max_subset_size,&dummy_idx,max_subset_size*max_subset_size,&work);CHKERRQ(ierr);
1059 
1060     /* matrices for adaptive selection */
1061     if (compute_Stilda) {
1062       if (!sub_schurs->sum_S_Ej_tilda_all) {
1063         ierr = MatCreate(PETSC_COMM_SELF,&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
1064         ierr = MatSetSizes(sub_schurs->sum_S_Ej_tilda_all,PETSC_DECIDE,PETSC_DECIDE,size_active_schur,size_active_schur);CHKERRQ(ierr);
1065         ierr = MatSetType(sub_schurs->sum_S_Ej_tilda_all,MATAIJ);CHKERRQ(ierr);
1066         ierr = MatSeqAIJSetPreallocation(sub_schurs->sum_S_Ej_tilda_all,0,nnz);CHKERRQ(ierr);
1067       }
1068       if (!sub_schurs->sum_S_Ej_inv_all) {
1069         ierr = MatCreate(PETSC_COMM_SELF,&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
1070         ierr = MatSetSizes(sub_schurs->sum_S_Ej_inv_all,PETSC_DECIDE,PETSC_DECIDE,size_active_schur,size_active_schur);CHKERRQ(ierr);
1071         ierr = MatSetType(sub_schurs->sum_S_Ej_inv_all,MATAIJ);CHKERRQ(ierr);
1072         ierr = MatSeqAIJSetPreallocation(sub_schurs->sum_S_Ej_inv_all,0,nnz);CHKERRQ(ierr);
1073       }
1074     }
1075 
1076     /* S_Ej_all */
1077     cum = cum2 = 0;
1078     ierr = MatDenseGetArray(S_all,&S_data);CHKERRQ(ierr);
1079     for (i=0;i<sub_schurs->n_subs;i++) {
1080       PetscInt j;
1081 
1082       /* get S_E */
1083       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
1084       if (S_lower_triangular) { /* I need to expand the upper triangular data (column oriented) */
1085         PetscInt k;
1086         for (k=0;k<subset_size;k++) {
1087           for (j=k;j<subset_size;j++) {
1088             work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
1089             work[j*subset_size+k] = PetscConj(S_data[cum2+k*size_schur+j]);
1090           }
1091         }
1092       } else { /* just copy to workspace */
1093         PetscInt k;
1094         for (k=0;k<subset_size;k++) {
1095           for (j=0;j<subset_size;j++) {
1096             work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
1097           }
1098         }
1099       }
1100       /* insert S_E values */
1101       for (j=0;j<subset_size;j++) dummy_idx[j] = cum+j;
1102       ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
1103 
1104       /* if adaptivity is requested, invert S_E blocks */
1105       if (compute_Stilda) {
1106         ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr);
1107         ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
1108         if (!use_getr) { /* TODO add sytrf/i for symmetric non hermitian */
1109           PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,work,&B_N,&B_ierr));
1110           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
1111           PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,work,&B_N,&B_ierr));
1112           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
1113         } else {
1114           PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,work,&B_N,pivots,&B_ierr));
1115           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
1116           PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,work,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
1117           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
1118         }
1119         ierr = PetscFPTrapPop();CHKERRQ(ierr);
1120         ierr = MatSetValues(sub_schurs->sum_S_Ej_inv_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
1121       }
1122       cum += subset_size;
1123       cum2 += subset_size*(size_schur + 1);
1124     }
1125     ierr = MatDenseRestoreArray(S_all,&S_data);CHKERRQ(ierr);
1126 
1127     if (solver_S) {
1128       ierr = MatFactorRestoreSchurComplement(F,&S_all);CHKERRQ(ierr);
1129     }
1130 
1131     schur_factor = NULL;
1132     if (compute_Stilda && size_active_schur) {
1133 
1134       if (sub_schurs->n_subs == 1 && size_schur == size_active_schur) { /* we already computed the inverse */
1135         PetscInt j;
1136         for (j=0;j<size_schur;j++) dummy_idx[j] = j;
1137         ierr = MatSetValues(sub_schurs->sum_S_Ej_tilda_all,size_schur,dummy_idx,size_schur,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
1138       } else {
1139         Mat S_all_inv=NULL;
1140         if (solver_S) { /* use MatFactor calls to invert S */
1141             /* let MatFactor factorize the Schur complement */
1142           ierr = MatFactorFactorizeSchurComplement(F);CHKERRQ(ierr);
1143           /* for adaptive selection we need S^-1; for solver reusage we need S_\Delta\Delta^-1.
1144              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 */
1145           if (factor_workaround) {
1146             PetscScalar *data;
1147             PetscInt nd = 0;
1148 
1149             if (!sub_schurs->is_posdef) {
1150               SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor update not yet implemented for non SPD matrices");
1151             }
1152             ierr = MatFactorGetSchurComplement(F,&S_all_inv);CHKERRQ(ierr);
1153             ierr = MatDenseGetArray(S_all_inv,&data);CHKERRQ(ierr);
1154             if (sub_schurs->is_dir) { /* dirichlet dofs could have different scalings */
1155               ierr = ISGetLocalSize(sub_schurs->is_dir,&nd);CHKERRQ(ierr);
1156             }
1157             ierr = PetscMalloc1((size_active_schur*(size_active_schur+1))/2+nd,&schur_factor);CHKERRQ(ierr);
1158             cum = 0;
1159             for (i=0;i<size_active_schur;i++) {
1160               ierr = PetscMemcpy(schur_factor+cum,data+i*(size_schur+1),(size_active_schur-i)*sizeof(PetscScalar));CHKERRQ(ierr);
1161               cum += size_active_schur-i;
1162             }
1163             for (i=0;i<nd;i++) schur_factor[cum+i] = data[(i+size_active_schur)*(size_schur+1)];
1164             /* invert without calling MatFactorInvertSchurComplement, since we are hacking */
1165             ierr = PetscBLASIntCast(size_schur,&B_N);CHKERRQ(ierr);
1166             ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
1167             PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,data,&B_N,&B_ierr));
1168             ierr = PetscFPTrapPop();CHKERRQ(ierr);
1169             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
1170             ierr = MatDenseRestoreArray(S_all_inv,&data);CHKERRQ(ierr);
1171           } else {
1172             ierr = MatFactorInvertSchurComplement(F);CHKERRQ(ierr);
1173             ierr = MatFactorGetSchurComplement(F,&S_all_inv);CHKERRQ(ierr);
1174           }
1175         } else { /* we need to invert explicitly since we are not using MatFactor for S */
1176           ierr = PetscObjectReference((PetscObject)S_all);CHKERRQ(ierr);
1177           S_all_inv = S_all;
1178           ierr = MatDenseGetArray(S_all_inv,&S_data);CHKERRQ(ierr);
1179           ierr = PetscBLASIntCast(size_schur,&B_N);CHKERRQ(ierr);
1180           ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
1181           if (!use_getr) { /* TODO add sytrf/i for symmetric non hermitian */
1182             PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,S_data,&B_N,&B_ierr));
1183             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
1184             PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,S_data,&B_N,&B_ierr));
1185             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
1186           } else {
1187             PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,S_data,&B_N,pivots,&B_ierr));
1188             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
1189             PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,S_data,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
1190             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
1191           }
1192           ierr = PetscFPTrapPop();CHKERRQ(ierr);
1193           ierr = MatDenseRestoreArray(S_all_inv,&S_data);CHKERRQ(ierr);
1194         }
1195         /* S_Ej_tilda_all */
1196         cum = cum2 = 0;
1197         ierr = MatDenseGetArray(S_all_inv,&S_data);CHKERRQ(ierr);
1198         for (i=0;i<sub_schurs->n_subs;i++) {
1199           PetscInt j;
1200 
1201           ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
1202           /* get (St^-1)_E */
1203           /* Here I don't need to expand to upper triangular since St^-1
1204              will be properly accessed later during adaptive selection */
1205           if (S_lower_triangular) {
1206             PetscInt k;
1207             for (k=0;k<subset_size;k++) {
1208               for (j=k;j<subset_size;j++) {
1209                 work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
1210               }
1211             }
1212           } else { /* copy to workspace */
1213             PetscInt k;
1214             for (k=0;k<subset_size;k++) {
1215               for (j=0;j<subset_size;j++) {
1216                 work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
1217               }
1218             }
1219           }
1220           for (j=0;j<subset_size;j++) dummy_idx[j] = cum+j;
1221           ierr = MatSetValues(sub_schurs->sum_S_Ej_tilda_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
1222           cum += subset_size;
1223           cum2 += subset_size*(size_schur + 1);
1224         }
1225         ierr = MatDenseRestoreArray(S_all_inv,&S_data);CHKERRQ(ierr);
1226         if (solver_S) {
1227           ierr = MatFactorRestoreSchurComplement(F,&S_all_inv);CHKERRQ(ierr);
1228         }
1229         ierr = MatDestroy(&S_all_inv);CHKERRQ(ierr);
1230       }
1231 
1232       /* move back factors to Schur data into F */
1233       if (factor_workaround) {
1234         Mat S_tmp;
1235         PetscInt nv,nd=0;
1236 
1237         if (!solver_S) {
1238           SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"This should not happen");
1239         }
1240         ierr = ISGetLocalSize(sub_schurs->is_vertices,&nv);CHKERRQ(ierr);
1241         ierr = MatFactorGetSchurComplement(F,&S_tmp);CHKERRQ(ierr);
1242         if (sub_schurs->is_posdef) {
1243           PetscScalar *data;
1244 
1245           ierr = MatZeroEntries(S_tmp);CHKERRQ(ierr);
1246           ierr = MatDenseGetArray(S_tmp,&data);CHKERRQ(ierr);
1247 
1248           if (S_lower_triangular) {
1249             cum = 0;
1250             for (i=0;i<size_active_schur;i++) {
1251               ierr = PetscMemcpy(data+i*(size_schur+1),schur_factor+cum,(size_active_schur-i)*sizeof(PetscScalar));CHKERRQ(ierr);
1252               cum += size_active_schur-i;
1253 	    }
1254           } else {
1255             ierr = PetscMemcpy(data,schur_factor,size_schur*size_schur*sizeof(PetscScalar));CHKERRQ(ierr);
1256           }
1257           if (sub_schurs->is_dir) {
1258             ierr = ISGetLocalSize(sub_schurs->is_dir,&nd);CHKERRQ(ierr);
1259 	    for (i=0;i<nd;i++) {
1260 	      data[(i+size_active_schur)*(size_schur+1)] = schur_factor[cum+i];
1261 	    }
1262           }
1263           /* workaround: since I cannot modify the matrices used inside the solvers for the forward and backward substitutions,
1264              set the diagonal entry of the Schur factor to a very large value */
1265           for (i=size_active_schur+nd;i<size_schur;i++) {
1266             data[i*(size_schur+1)] = PETSC_MAX_REAL;
1267           }
1268           ierr = MatDenseRestoreArray(S_tmp,&data);CHKERRQ(ierr);
1269         } else {
1270           SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor update not yet implemented for non SPD matrices");
1271         }
1272         ierr = MatFactorRestoreSchurComplement(F,&S_tmp);CHKERRQ(ierr);
1273       }
1274     }
1275     ierr = PetscFree2(dummy_idx,work);CHKERRQ(ierr);
1276     ierr = PetscFree(schur_factor);CHKERRQ(ierr);
1277   }
1278   ierr = PetscFree(nnz);CHKERRQ(ierr);
1279   ierr = MatDestroy(&F);CHKERRQ(ierr);
1280   ierr = ISDestroy(&is_I_layer);CHKERRQ(ierr);
1281   ierr = MatDestroy(&S_all);CHKERRQ(ierr);
1282   ierr = MatDestroy(&A_BB);CHKERRQ(ierr);
1283   ierr = MatDestroy(&A_IB);CHKERRQ(ierr);
1284   ierr = MatDestroy(&A_BI);CHKERRQ(ierr);
1285   ierr = MatAssemblyBegin(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1286   ierr = MatAssemblyEnd(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1287   if (compute_Stilda) {
1288     ierr = MatAssemblyBegin(sub_schurs->sum_S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1289     ierr = MatAssemblyEnd(sub_schurs->sum_S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1290     ierr = MatAssemblyBegin(sub_schurs->sum_S_Ej_inv_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1291     ierr = MatAssemblyEnd(sub_schurs->sum_S_Ej_inv_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1292   }
1293 
1294   /* Global matrix of all assembled Schur on subsets */
1295   ierr = MatISSetLocalMat(work_mat,sub_schurs->S_Ej_all);CHKERRQ(ierr);
1296   ierr = MatISSetMPIXAIJPreallocation_Private(work_mat,global_schur_subsets,PETSC_TRUE);CHKERRQ(ierr);
1297   ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
1298 
1299   /* Get local part of (\sum_j S_Ej) */
1300   if (!sub_schurs->sum_S_Ej_all) {
1301     ierr = MatGetSubMatrices(global_schur_subsets,1,&all_subsets_n,&all_subsets_n,MAT_INITIAL_MATRIX,&submats);CHKERRQ(ierr);
1302     sub_schurs->sum_S_Ej_all = submats[0];
1303   } else {
1304     ierr = PetscMalloc1(1,&submats);CHKERRQ(ierr);
1305     submats[0] = sub_schurs->sum_S_Ej_all;
1306     ierr = MatGetSubMatrices(global_schur_subsets,1,&all_subsets_n,&all_subsets_n,MAT_REUSE_MATRIX,&submats);CHKERRQ(ierr);
1307   }
1308 
1309   /* Compute explicitly (\sum_j S_Ej)^-1 (faster scaling during PCApply, needs extra work when doing setup) */
1310   if (faster_deluxe) {
1311     Mat         tmpmat;
1312     PetscScalar *array;
1313     PetscInt    cum;
1314 
1315     ierr = MatSeqAIJGetArray(sub_schurs->sum_S_Ej_all,&array);CHKERRQ(ierr);
1316     cum = 0;
1317     for (i=0;i<sub_schurs->n_subs;i++) {
1318       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
1319       ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr);
1320       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
1321       if (!use_getr) {
1322         PetscInt j,k;
1323 
1324         PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,array+cum,&B_N,&B_ierr));
1325         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
1326         PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,array+cum,&B_N,&B_ierr));
1327         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
1328         for (j=0;j<B_N;j++) {
1329           for (k=j+1;k<B_N;k++) {
1330             array[k*B_N+j+cum] = array[j*B_N+k+cum];
1331           }
1332         }
1333       } else {
1334         PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,array+cum,&B_N,pivots,&B_ierr));
1335         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
1336         PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,array+cum,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
1337         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
1338       }
1339       ierr = PetscFPTrapPop();CHKERRQ(ierr);
1340       cum += subset_size*subset_size;
1341     }
1342     ierr = MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_all,&array);CHKERRQ(ierr);
1343     ierr = MatMatMult(sub_schurs->S_Ej_all,sub_schurs->sum_S_Ej_all,MAT_INITIAL_MATRIX,1.0,&tmpmat);CHKERRQ(ierr);
1344     ierr = MatDestroy(&sub_schurs->S_Ej_all);CHKERRQ(ierr);
1345     ierr = MatDestroy(&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
1346     sub_schurs->S_Ej_all = tmpmat;
1347   }
1348 
1349   /* Get local part of (\sum_j S^-1_Ej) (\sum_j St^-1_Ej) */
1350   if (compute_Stilda) {
1351     ierr = MatISSetLocalMat(work_mat,sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
1352     ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
1353     submats[0] = sub_schurs->sum_S_Ej_tilda_all;
1354     ierr = MatGetSubMatrices(global_schur_subsets,1,&all_subsets_n,&all_subsets_n,MAT_REUSE_MATRIX,&submats);CHKERRQ(ierr);
1355     ierr = MatISSetLocalMat(work_mat,sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
1356     ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
1357     submats[0] = sub_schurs->sum_S_Ej_inv_all;
1358     ierr = MatGetSubMatrices(global_schur_subsets,1,&all_subsets_n,&all_subsets_n,MAT_REUSE_MATRIX,&submats);CHKERRQ(ierr);
1359   }
1360 
1361   /* free workspace */
1362   ierr = PetscFree(submats);CHKERRQ(ierr);
1363   ierr = PetscFree2(Bwork,pivots);CHKERRQ(ierr);
1364   ierr = MatDestroy(&global_schur_subsets);CHKERRQ(ierr);
1365   ierr = MatDestroy(&work_mat);CHKERRQ(ierr);
1366   ierr = ISDestroy(&all_subsets_n);CHKERRQ(ierr);
1367   ierr = PetscCommDestroy(&comm_n);CHKERRQ(ierr);
1368   PetscFunctionReturn(0);
1369 }
1370 
1371 #undef __FUNCT__
1372 #define __FUNCT__ "PCBDDCSubSchursInit"
1373 PetscErrorCode PCBDDCSubSchursInit(PCBDDCSubSchurs sub_schurs, IS is_I, IS is_B, PCBDDCGraph graph, ISLocalToGlobalMapping BtoNmap)
1374 {
1375   IS              *faces,*edges,*all_cc,vertices;
1376   PetscInt        i,n_faces,n_edges,n_all_cc;
1377   PetscBool       is_sorted;
1378   PetscErrorCode  ierr;
1379 
1380   PetscFunctionBegin;
1381   ierr = ISSorted(is_I,&is_sorted);CHKERRQ(ierr);
1382   if (!is_sorted) {
1383     SETERRQ(PetscObjectComm((PetscObject)is_I),PETSC_ERR_PLIB,"IS for I dofs should be shorted");
1384   }
1385   ierr = ISSorted(is_B,&is_sorted);CHKERRQ(ierr);
1386   if (!is_sorted) {
1387     SETERRQ(PetscObjectComm((PetscObject)is_B),PETSC_ERR_PLIB,"IS for B dofs should be shorted");
1388   }
1389 
1390   /* reset any previous data */
1391   ierr = PCBDDCSubSchursReset(sub_schurs);CHKERRQ(ierr);
1392 
1393   /* get index sets for faces and edges (already sorted by global ordering) */
1394   ierr = PCBDDCGraphGetCandidatesIS(graph,&n_faces,&faces,&n_edges,&edges,&vertices);CHKERRQ(ierr);
1395   n_all_cc = n_faces+n_edges;
1396   ierr = PetscBTCreate(n_all_cc,&sub_schurs->is_edge);CHKERRQ(ierr);
1397   ierr = PetscMalloc1(n_all_cc,&all_cc);CHKERRQ(ierr);
1398   for (i=0;i<n_faces;i++) {
1399     all_cc[i] = faces[i];
1400   }
1401   for (i=0;i<n_edges;i++) {
1402     all_cc[n_faces+i] = edges[i];
1403     ierr = PetscBTSet(sub_schurs->is_edge,n_faces+i);CHKERRQ(ierr);
1404   }
1405   ierr = PetscFree(faces);CHKERRQ(ierr);
1406   ierr = PetscFree(edges);CHKERRQ(ierr);
1407   sub_schurs->is_dir = NULL;
1408   ierr = PCBDDCGraphGetDirichletDofsB(graph,&sub_schurs->is_dir);CHKERRQ(ierr);
1409 
1410   /* Determine if MatFactor can be used */
1411   sub_schurs->schur_explicit = PETSC_FALSE;
1412 #if defined(PETSC_HAVE_MUMPS)
1413   sub_schurs->schur_explicit = PETSC_TRUE;
1414 #endif
1415 #if defined(PETSC_HAVE_MKL_PARDISO)
1416   sub_schurs->schur_explicit = PETSC_TRUE;
1417 #endif
1418 
1419   ierr = PetscObjectReference((PetscObject)is_I);CHKERRQ(ierr);
1420   sub_schurs->is_I = is_I;
1421   ierr = PetscObjectReference((PetscObject)is_B);CHKERRQ(ierr);
1422   sub_schurs->is_B = is_B;
1423   ierr = PetscObjectReference((PetscObject)graph->l2gmap);CHKERRQ(ierr);
1424   sub_schurs->l2gmap = graph->l2gmap;
1425   ierr = PetscObjectReference((PetscObject)BtoNmap);CHKERRQ(ierr);
1426   sub_schurs->BtoNmap = BtoNmap;
1427   sub_schurs->n_subs = n_all_cc;
1428   sub_schurs->is_subs = all_cc;
1429   if (!sub_schurs->schur_explicit) { /* sort by local ordering if we are not using MatFactor */
1430     for (i=0;i<sub_schurs->n_subs;i++) {
1431       ierr = ISSort(sub_schurs->is_subs[i]);CHKERRQ(ierr);
1432     }
1433   }
1434   sub_schurs->is_vertices = vertices;
1435   sub_schurs->S_Ej_all = NULL;
1436   sub_schurs->sum_S_Ej_all = NULL;
1437   sub_schurs->sum_S_Ej_inv_all = NULL;
1438   sub_schurs->sum_S_Ej_tilda_all = NULL;
1439   sub_schurs->is_Ej_all = NULL;
1440   PetscFunctionReturn(0);
1441 }
1442 
1443 #undef __FUNCT__
1444 #define __FUNCT__ "PCBDDCSubSchursCreate"
1445 PetscErrorCode PCBDDCSubSchursCreate(PCBDDCSubSchurs *sub_schurs)
1446 {
1447   PCBDDCSubSchurs schurs_ctx;
1448   PetscErrorCode  ierr;
1449 
1450   PetscFunctionBegin;
1451   ierr = PetscNew(&schurs_ctx);CHKERRQ(ierr);
1452   schurs_ctx->n_subs = 0;
1453   *sub_schurs = schurs_ctx;
1454   PetscFunctionReturn(0);
1455 }
1456 
1457 #undef __FUNCT__
1458 #define __FUNCT__ "PCBDDCSubSchursDestroy"
1459 PetscErrorCode PCBDDCSubSchursDestroy(PCBDDCSubSchurs *sub_schurs)
1460 {
1461   PetscErrorCode ierr;
1462 
1463   PetscFunctionBegin;
1464   ierr = PCBDDCSubSchursReset(*sub_schurs);CHKERRQ(ierr);
1465   ierr = PetscFree(*sub_schurs);CHKERRQ(ierr);
1466   PetscFunctionReturn(0);
1467 }
1468 
1469 #undef __FUNCT__
1470 #define __FUNCT__ "PCBDDCSubSchursReset"
1471 PetscErrorCode PCBDDCSubSchursReset(PCBDDCSubSchurs sub_schurs)
1472 {
1473   PetscInt       i;
1474   PetscErrorCode ierr;
1475 
1476   PetscFunctionBegin;
1477   ierr = MatDestroy(&sub_schurs->A);CHKERRQ(ierr);
1478   ierr = MatDestroy(&sub_schurs->S);CHKERRQ(ierr);
1479   ierr = ISDestroy(&sub_schurs->is_I);CHKERRQ(ierr);
1480   ierr = ISDestroy(&sub_schurs->is_B);CHKERRQ(ierr);
1481   ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->l2gmap);CHKERRQ(ierr);
1482   ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->BtoNmap);CHKERRQ(ierr);
1483   ierr = MatDestroy(&sub_schurs->S_Ej_all);CHKERRQ(ierr);
1484   ierr = MatDestroy(&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
1485   ierr = MatDestroy(&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
1486   ierr = MatDestroy(&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
1487   ierr = ISDestroy(&sub_schurs->is_Ej_all);CHKERRQ(ierr);
1488   ierr = ISDestroy(&sub_schurs->is_vertices);CHKERRQ(ierr);
1489   ierr = ISDestroy(&sub_schurs->is_dir);CHKERRQ(ierr);
1490   ierr = PetscBTDestroy(&sub_schurs->is_edge);CHKERRQ(ierr);
1491   for (i=0;i<sub_schurs->n_subs;i++) {
1492     ierr = ISDestroy(&sub_schurs->is_subs[i]);CHKERRQ(ierr);
1493   }
1494   if (sub_schurs->n_subs) {
1495     ierr = PetscFree(sub_schurs->is_subs);CHKERRQ(ierr);
1496   }
1497   if (sub_schurs->reuse_solver) {
1498     ierr = PCBDDCReuseSolversReset(sub_schurs->reuse_solver);CHKERRQ(ierr);
1499   }
1500   ierr = PetscFree(sub_schurs->reuse_solver);CHKERRQ(ierr);
1501   sub_schurs->n_subs = 0;
1502   PetscFunctionReturn(0);
1503 }
1504 
1505 #undef __FUNCT__
1506 #define __FUNCT__ "PCBDDCAdjGetNextLayer_Private"
1507 PETSC_STATIC_INLINE PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt* queue_tip,PetscInt n_prev,PetscBT touched,PetscInt* xadj,PetscInt* adjncy,PetscInt* n_added)
1508 {
1509   PetscInt       i,j,n;
1510   PetscErrorCode ierr;
1511 
1512   PetscFunctionBegin;
1513   n = 0;
1514   for (i=-n_prev;i<0;i++) {
1515     PetscInt start_dof = queue_tip[i];
1516     for (j=xadj[start_dof];j<xadj[start_dof+1];j++) {
1517       PetscInt dof = adjncy[j];
1518       if (!PetscBTLookup(touched,dof)) {
1519         ierr = PetscBTSet(touched,dof);CHKERRQ(ierr);
1520         queue_tip[n] = dof;
1521         n++;
1522       }
1523     }
1524   }
1525   *n_added = n;
1526   PetscFunctionReturn(0);
1527 }
1528