xref: /petsc/src/ksp/pc/impls/bddc/bddcschurs.c (revision 5ec10c6ac20a3057719a2d68a0b1ec85fc87085b)
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 
8 #undef __FUNCT__
9 #define __FUNCT__ "PCBDDCComputeExplicitSchur"
10 static PetscErrorCode PCBDDCComputeExplicitSchur(Mat M, PetscBool issym, MatReuse reuse, Mat *S)
11 {
12   Mat            B, C, D, Bd, Cd, AinvBd;
13   KSP            ksp;
14   PC             pc;
15   PetscBool      isLU, isILU, isCHOL, Bdense, Cdense;
16   PetscReal      fill = 2.0;
17   PetscMPIInt    size;
18   PetscErrorCode ierr;
19 
20   PetscFunctionBegin;
21   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)M),&size);CHKERRQ(ierr);
22   if (size != 1) {
23     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not for parallel matrices");
24   }
25   ierr = MatSchurComplementGetSubMatrices(M, NULL, NULL, &B, &C, &D);CHKERRQ(ierr);
26   ierr = MatSchurComplementGetKSP(M, &ksp);CHKERRQ(ierr);
27   ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr);
28   ierr = PetscObjectTypeCompare((PetscObject) pc, PCLU, &isLU);CHKERRQ(ierr);
29   ierr = PetscObjectTypeCompare((PetscObject) pc, PCILU, &isILU);CHKERRQ(ierr);
30   ierr = PetscObjectTypeCompare((PetscObject) pc, PCCHOLESKY, &isCHOL);CHKERRQ(ierr);
31   ierr = PetscObjectTypeCompare((PetscObject) B, MATSEQDENSE, &Bdense);CHKERRQ(ierr);
32   ierr = PetscObjectTypeCompare((PetscObject) C, MATSEQDENSE, &Cdense);CHKERRQ(ierr);
33   if (!Bdense) {
34     ierr = MatConvert(B, MATSEQDENSE, MAT_INITIAL_MATRIX, &Bd);CHKERRQ(ierr);
35   } else {
36     Bd = B;
37   }
38 
39   if (isLU || isILU || isCHOL) {
40     Mat fact;
41 
42     ierr = KSPSetUp(ksp);CHKERRQ(ierr);
43     ierr = PCFactorGetMatrix(pc, &fact);CHKERRQ(ierr);
44     ierr = MatDuplicate(Bd, MAT_DO_NOT_COPY_VALUES, &AinvBd);CHKERRQ(ierr);
45     ierr = MatMatSolve(fact, Bd, AinvBd);CHKERRQ(ierr);
46   } else {
47     Mat Ainvd;
48 
49     ierr = PCComputeExplicitOperator(pc, &Ainvd);CHKERRQ(ierr);
50     ierr = MatMatMult(Ainvd, Bd, MAT_INITIAL_MATRIX, fill, &AinvBd);CHKERRQ(ierr);
51     ierr = MatDestroy(&Ainvd);CHKERRQ(ierr);
52   }
53   if (!Bdense & !issym) {
54     ierr = MatDestroy(&Bd);CHKERRQ(ierr);
55   }
56 
57   if (!issym) {
58     if (!Cdense) {
59       ierr = MatConvert(C, MATSEQDENSE, MAT_INITIAL_MATRIX, &Cd);CHKERRQ(ierr);
60     } else {
61       Cd = C;
62     }
63     ierr = MatMatMult(Cd, AinvBd, reuse, fill, S);CHKERRQ(ierr);
64     if (!Cdense) {
65       ierr = MatDestroy(&Cd);CHKERRQ(ierr);
66     }
67   } else {
68     ierr = MatTransposeMatMult(Bd, AinvBd, reuse, fill, S);CHKERRQ(ierr);
69     if (!Bdense) {
70       ierr = MatDestroy(&Bd);CHKERRQ(ierr);
71     }
72   }
73   ierr = MatDestroy(&AinvBd);CHKERRQ(ierr);
74 
75   if (D) {
76     Mat       Dd;
77     PetscBool Ddense;
78 
79     ierr = PetscObjectTypeCompare((PetscObject)D,MATSEQDENSE,&Ddense);CHKERRQ(ierr);
80     if (!Ddense) {
81       ierr = MatConvert(D, MATSEQDENSE, MAT_INITIAL_MATRIX, &Dd);CHKERRQ(ierr);
82     } else {
83       Dd = D;
84     }
85     ierr = MatAYPX(*S,-1.0,Dd,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
86     if (!Ddense) {
87       ierr = MatDestroy(&Dd);CHKERRQ(ierr);
88     }
89   } else {
90     ierr = MatScale(*S,-1.0);CHKERRQ(ierr);
91   }
92   PetscFunctionReturn(0);
93 }
94 
95 #undef __FUNCT__
96 #define __FUNCT__ "PCBDDCSubSchursSetUp"
97 PetscErrorCode PCBDDCSubSchursSetUp(PCBDDCSubSchurs sub_schurs, PetscInt xadj[], PetscInt adjncy[], PetscInt nlayers, PetscBool compute_Stilda, PetscBool deluxe, PetscBool invert_Stildas, PetscBool use_edges, PetscBool use_faces)
98 {
99   Mat                    A_II,A_IB,A_BI,A_BB;
100   Mat                    AE_II,*AE_IE,*AE_EI,*AE_EE;
101   Mat                    S_all,global_schur_subsets,work_mat;
102   Mat                    S_Ej_tilda_all,S_Ej_inv_all;
103   ISLocalToGlobalMapping l2gmap_subsets;
104   IS                     is_I,*is_subset_B,temp_is;
105   PetscInt               *nnz,*all_local_idx_G,*all_local_idx_B,*all_local_idx_N;
106   PetscInt               i,subset_size,max_subset_size_seq;
107   PetscInt               extra,local_size,global_size;
108   PetscBLASInt           B_N,B_ierr,B_lwork,*pivots;
109   PetscScalar            *work;
110   PetscErrorCode         ierr;
111 
112   PetscFunctionBegin;
113   /* get Schur complement matrices */
114   if (!sub_schurs->use_mumps) {
115     if (compute_Stilda) {
116       SETERRQ(PetscObjectComm((PetscObject)sub_schurs->l2gmap),PETSC_ERR_SUP,"Adaptive selection of constraints requires MUMPS");
117     }
118     ierr = MatSchurComplementGetSubMatrices(sub_schurs->S,&A_II,NULL,&A_IB,&A_BI,&A_BB);CHKERRQ(ierr);
119     ierr = PetscMalloc4(sub_schurs->n_subs,&is_subset_B,
120                         sub_schurs->n_subs,&AE_IE,
121                         sub_schurs->n_subs,&AE_EI,
122                         sub_schurs->n_subs,&AE_EE);CHKERRQ(ierr);
123   } else {
124     is_subset_B = NULL;
125     AE_IE = NULL;
126     AE_EI = NULL;
127     AE_EE = NULL;
128   }
129 
130   /* determine interior problems */
131   ierr = ISDestroy(&sub_schurs->is_I_layer);CHKERRQ(ierr);
132   ierr = ISGetLocalSize(sub_schurs->is_I,&i);CHKERRQ(ierr);
133   if (nlayers >= 0 && i) { /* Interior problems can be different from the original one */
134     PetscBT                touched;
135     const PetscInt*        idx_B;
136     PetscInt               n_I,n_B,n_local_dofs,n_prev_added,j,layer,*local_numbering;
137 
138     if (xadj == NULL || adjncy == NULL) {
139       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cannot request layering without adjacency");
140     }
141     /* get sizes */
142     ierr = ISGetLocalSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr);
143     ierr = ISGetLocalSize(sub_schurs->is_B,&n_B);CHKERRQ(ierr);
144 
145     ierr = PetscMalloc1(n_I+n_B,&local_numbering);CHKERRQ(ierr);
146     ierr = PetscBTCreate(n_I+n_B,&touched);CHKERRQ(ierr);
147     ierr = PetscBTMemzero(n_I+n_B,touched);CHKERRQ(ierr);
148 
149     /* all boundary dofs must be skipped when adding layers */
150     ierr = ISGetIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr);
151     for (j=0;j<n_B;j++) {
152       ierr = PetscBTSet(touched,idx_B[j]);CHKERRQ(ierr);
153     }
154     ierr = PetscMemcpy(local_numbering,idx_B,n_B*sizeof(PetscInt));CHKERRQ(ierr);
155     ierr = ISRestoreIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr);
156 
157     /* add prescribed number of layers of dofs */
158     n_local_dofs = n_B;
159     n_prev_added = n_B;
160     for (layer=0;layer<nlayers;layer++) {
161       PetscInt n_added;
162       if (n_local_dofs == n_I+n_B) break;
163       if (n_local_dofs > n_I+n_B) {
164         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);
165       }
166       ierr = PCBDDCAdjGetNextLayer_Private(local_numbering+n_local_dofs,n_prev_added,touched,xadj,adjncy,&n_added);CHKERRQ(ierr);
167       n_prev_added = n_added;
168       n_local_dofs += n_added;
169       if (!n_added) break;
170     }
171     ierr = PetscBTDestroy(&touched);CHKERRQ(ierr);
172 
173     /* IS for I layer dofs in original numbering */
174     ierr = ISCreateGeneral(PetscObjectComm((PetscObject)sub_schurs->is_I),n_local_dofs-n_B,local_numbering+n_B,PETSC_COPY_VALUES,&sub_schurs->is_I_layer);CHKERRQ(ierr);
175     ierr = PetscFree(local_numbering);CHKERRQ(ierr);
176     ierr = ISSort(sub_schurs->is_I_layer);CHKERRQ(ierr);
177     /* IS for I layer dofs in I numbering */
178     if (!sub_schurs->use_mumps) {
179       ISLocalToGlobalMapping ItoNmap;
180       ierr = ISLocalToGlobalMappingCreateIS(sub_schurs->is_I,&ItoNmap);CHKERRQ(ierr);
181       ierr = ISGlobalToLocalMappingApplyIS(ItoNmap,IS_GTOLM_DROP,sub_schurs->is_I_layer,&is_I);CHKERRQ(ierr);
182       ierr = ISLocalToGlobalMappingDestroy(&ItoNmap);CHKERRQ(ierr);
183 
184       /* II block */
185       ierr = MatGetSubMatrix(A_II,is_I,is_I,MAT_INITIAL_MATRIX,&AE_II);CHKERRQ(ierr);
186     }
187   } else {
188     PetscInt n_I;
189 
190     /* IS for I dofs in original numbering */
191     ierr = PetscObjectReference((PetscObject)sub_schurs->is_I);CHKERRQ(ierr);
192     sub_schurs->is_I_layer = sub_schurs->is_I;
193 
194     /* IS for I dofs in I numbering (strided 1) */
195     if (!sub_schurs->use_mumps) {
196       ierr = ISGetSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr);
197       ierr = ISCreateStride(PetscObjectComm((PetscObject)sub_schurs->is_I),n_I,0,1,&is_I);CHKERRQ(ierr);
198 
199       /* II block is the same */
200       ierr = PetscObjectReference((PetscObject)A_II);CHKERRQ(ierr);
201       AE_II = A_II;
202     }
203   }
204   /* Get info on subset sizes and sum of all subsets sizes */
205   max_subset_size_seq = 0;
206   local_size = 0;
207   for (i=0;i<sub_schurs->n_subs_seq;i++) {
208     PetscInt j = sub_schurs->index_sequential[i];
209     ierr = ISGetLocalSize(sub_schurs->is_subs[j],&subset_size);CHKERRQ(ierr);
210     max_subset_size_seq = PetscMax(subset_size,max_subset_size_seq);
211     local_size += subset_size;
212   }
213 
214   /* Work arrays for local indices */
215   ierr = PetscMalloc1(local_size,&all_local_idx_B);CHKERRQ(ierr);
216   extra = 0;
217   if (sub_schurs->use_mumps) {
218     ierr = ISGetLocalSize(sub_schurs->is_I_layer,&extra);CHKERRQ(ierr);
219   }
220   ierr = PetscMalloc1(local_size+extra,&all_local_idx_N);CHKERRQ(ierr);
221   if (extra) {
222     const PetscInt *idxs;
223     ierr = ISGetIndices(sub_schurs->is_I_layer,&idxs);CHKERRQ(ierr);
224     ierr = PetscMemcpy(all_local_idx_N,idxs,extra*sizeof(PetscInt));CHKERRQ(ierr);
225     ierr = ISRestoreIndices(sub_schurs->is_I_layer,&idxs);CHKERRQ(ierr);
226   }
227   ierr = PetscMalloc1(local_size,&nnz);CHKERRQ(ierr);
228 
229   /* Get local indices in local numbering */
230   local_size = 0;
231   for (i=0;i<sub_schurs->n_subs_seq;i++) {
232     PetscInt j;
233     const    PetscInt *idxs;
234 
235     PetscInt local_problem_index = sub_schurs->index_sequential[i];
236     ierr = ISGetLocalSize(sub_schurs->is_subs[local_problem_index],&subset_size);CHKERRQ(ierr);
237     ierr = ISGetIndices(sub_schurs->is_subs[local_problem_index],&idxs);CHKERRQ(ierr);
238     /* subset indices in local numbering */
239     ierr = PetscMemcpy(all_local_idx_N+local_size+extra,idxs,subset_size*sizeof(PetscInt));CHKERRQ(ierr);
240     ierr = ISRestoreIndices(sub_schurs->is_subs[local_problem_index],&idxs);CHKERRQ(ierr);
241     for (j=0;j<subset_size;j++) nnz[local_size+j] = subset_size;
242     local_size += subset_size;
243   }
244 
245   S_all = NULL;
246   sub_schurs->is_hermitian = PETSC_FALSE;
247   sub_schurs->is_posdef = PETSC_FALSE;
248   if (sub_schurs->A) {
249     ierr = MatIsHermitian(sub_schurs->A,0.0,&sub_schurs->is_hermitian);CHKERRQ(ierr);
250     if (sub_schurs->is_hermitian) {
251       PetscScalar val;
252       Vec         vec1,vec2;
253 
254       ierr = MatCreateVecs(sub_schurs->A,&vec1,&vec2);CHKERRQ(ierr);
255       ierr = VecSetRandom(vec1,NULL);
256       ierr = VecCopy(vec1,vec2);CHKERRQ(ierr);
257       ierr = MatMult(sub_schurs->A,vec2,vec1);CHKERRQ(ierr);
258       ierr = VecDot(vec1,vec2,&val);CHKERRQ(ierr);
259       if (PetscRealPart(val) > 0. && PetscImaginaryPart(val) == 0.) sub_schurs->is_posdef = PETSC_TRUE;
260       ierr = VecDestroy(&vec1);CHKERRQ(ierr);
261       ierr = VecDestroy(&vec2);CHKERRQ(ierr);
262     }
263   }
264   if (sub_schurs->n_subs) {
265     if (!sub_schurs->use_mumps) {
266       /* subsets in original and boundary numbering */
267       for (i=0;i<sub_schurs->n_subs;i++) {
268         ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,sub_schurs->is_subs[i],&is_subset_B[i]);CHKERRQ(ierr);
269       }
270 
271       /* EE block */
272       for (i=0;i<sub_schurs->n_subs;i++) {
273         ierr = MatGetSubMatrix(A_BB,is_subset_B[i],is_subset_B[i],MAT_INITIAL_MATRIX,&AE_EE[i]);CHKERRQ(ierr);
274       }
275       /* IE block */
276       for (i=0;i<sub_schurs->n_subs;i++) {
277         ierr = MatGetSubMatrix(A_IB,is_I,is_subset_B[i],MAT_INITIAL_MATRIX,&AE_IE[i]);CHKERRQ(ierr);
278       }
279       /* EI block */
280       for (i=0;i<sub_schurs->n_subs;i++) {
281         ierr = MatGetSubMatrix(A_BI,is_subset_B[i],is_I,MAT_INITIAL_MATRIX,&AE_EI[i]);CHKERRQ(ierr);
282       }
283 
284       /* setup Schur complements on subset */
285       for (i=0;i<sub_schurs->n_subs;i++) {
286         ierr = MatDestroy(&sub_schurs->S_Ej[i]);CHKERRQ(ierr);
287         ierr = MatCreateSchurComplement(AE_II,AE_II,AE_IE[i],AE_EI[i],AE_EE[i],&sub_schurs->S_Ej[i]);CHKERRQ(ierr);
288         ierr = MatDestroy(&AE_EE[i]);CHKERRQ(ierr);
289         ierr = MatDestroy(&AE_IE[i]);CHKERRQ(ierr);
290         ierr = MatDestroy(&AE_EI[i]);CHKERRQ(ierr);
291         if (AE_II == A_II) { /* we can reuse the same ksp */
292           KSP ksp;
293           ierr = MatSchurComplementGetKSP(sub_schurs->S,&ksp);CHKERRQ(ierr);
294           ierr = MatSchurComplementSetKSP(sub_schurs->S_Ej[i],ksp);CHKERRQ(ierr);
295         } else { /* build new ksp object which inherits ksp and pc types from the original one */
296           KSP      origksp,schurksp;
297           PC       origpc,schurpc;
298           KSPType  ksp_type;
299           PCType   pc_type;
300           PetscInt n_internal;
301 
302           ierr = MatSchurComplementGetKSP(sub_schurs->S,&origksp);CHKERRQ(ierr);
303           ierr = MatSchurComplementGetKSP(sub_schurs->S_Ej[i],&schurksp);CHKERRQ(ierr);
304           ierr = KSPGetType(origksp,&ksp_type);CHKERRQ(ierr);
305           ierr = KSPSetType(schurksp,ksp_type);CHKERRQ(ierr);
306           ierr = KSPGetPC(schurksp,&schurpc);CHKERRQ(ierr);
307           ierr = KSPGetPC(origksp,&origpc);CHKERRQ(ierr);
308           ierr = PCGetType(origpc,&pc_type);CHKERRQ(ierr);
309           ierr = PCSetType(schurpc,pc_type);CHKERRQ(ierr);
310           ierr = ISGetSize(is_I,&n_internal);CHKERRQ(ierr);
311           if (n_internal) { /* UMFPACK gives error with 0 sized problems */
312             MatSolverPackage solver=NULL;
313             ierr = PCFactorGetMatSolverPackage(origpc,(const MatSolverPackage*)&solver);CHKERRQ(ierr);
314             if (solver) {
315               ierr = PCFactorSetMatSolverPackage(schurpc,solver);CHKERRQ(ierr);
316             }
317           }
318           ierr = KSPSetUp(schurksp);CHKERRQ(ierr);
319         }
320       }
321       /* free */
322       ierr = ISDestroy(&is_I);CHKERRQ(ierr);
323       ierr = MatDestroy(&AE_II);CHKERRQ(ierr);
324       for (i=0;i<sub_schurs->n_subs;i++) {
325         ierr = ISDestroy(&is_subset_B[i]);CHKERRQ(ierr);
326       }
327       ierr = PetscFree4(is_subset_B,AE_IE,AE_EI,AE_EE);CHKERRQ(ierr);
328 #if defined(PETSC_HAVE_MUMPS)
329     } else {
330       Mat           A,F;
331       IS            is_A_all;
332       PetscInt      *idxs_schur,n_I;
333 
334       /* get working mat */
335       ierr = ISGetLocalSize(sub_schurs->is_I_layer,&n_I);CHKERRQ(ierr);
336       ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size+n_I,all_local_idx_N,PETSC_COPY_VALUES,&is_A_all);CHKERRQ(ierr);
337       ierr = MatGetSubMatrixUnsorted(sub_schurs->A,is_A_all,is_A_all,MAT_INITIAL_MATRIX,&A);CHKERRQ(ierr);
338       ierr = ISDestroy(&is_A_all);CHKERRQ(ierr);
339 
340       if (n_I) {
341         if (sub_schurs->is_hermitian && sub_schurs->is_posdef) {
342           ierr = MatGetFactor(A,MATSOLVERMUMPS,MAT_FACTOR_CHOLESKY,&F);CHKERRQ(ierr);
343         } else {
344           ierr = MatGetFactor(A,MATSOLVERMUMPS,MAT_FACTOR_LU,&F);CHKERRQ(ierr);
345         }
346 
347         /* subsets ordered last */
348         ierr = PetscMalloc1(local_size,&idxs_schur);CHKERRQ(ierr);
349         for (i=0;i<local_size;i++) {
350           idxs_schur[i] = n_I+i+1;
351         }
352         ierr = MatMumpsSetSchurIndices(F,local_size,idxs_schur);CHKERRQ(ierr);
353         ierr = PetscFree(idxs_schur);CHKERRQ(ierr);
354 
355         /* factorization step */
356         if (sub_schurs->is_hermitian && sub_schurs->is_posdef) {
357           ierr = MatCholeskyFactorSymbolic(F,A,NULL,NULL);CHKERRQ(ierr);
358           ierr = MatCholeskyFactorNumeric(F,A,NULL);CHKERRQ(ierr);
359         } else {
360           ierr = MatLUFactorSymbolic(F,A,NULL,NULL,NULL);CHKERRQ(ierr);
361           ierr = MatLUFactorNumeric(F,A,NULL);CHKERRQ(ierr);
362         }
363 
364         /* get explicit Schur Complement computed during numeric factorization */
365         ierr = MatMumpsGetSchurComplement(F,&S_all);CHKERRQ(ierr);
366         ierr = MatDestroy(&F);CHKERRQ(ierr);
367       } else {
368         ierr = MatConvert(A,MATSEQDENSE,MAT_INITIAL_MATRIX,&S_all);CHKERRQ(ierr);
369       }
370       ierr = MatDestroy(&A);CHKERRQ(ierr);
371 #endif
372     }
373   } else {
374     ierr = PetscFree(nnz);CHKERRQ(ierr);
375     ierr = PetscFree(all_local_idx_N);CHKERRQ(ierr);
376   }
377 
378   if (!sub_schurs->n_subs_seq_g) {
379     PetscFunctionReturn(0);
380   }
381 
382   /* subset indices in local boundary numbering */
383   if (!sub_schurs->is_Ej_all) {
384     ierr = ISGlobalToLocalMappingApply(sub_schurs->BtoNmap,IS_GTOLM_DROP,local_size,all_local_idx_N+extra,&subset_size,all_local_idx_B);CHKERRQ(ierr);
385     if (subset_size != local_size) {
386       SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in sub_schurs serial (BtoNmap)! %d != %d\n",subset_size,local_size);
387     }
388     ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size,all_local_idx_B,PETSC_OWN_POINTER,&sub_schurs->is_Ej_all);CHKERRQ(ierr);
389   }
390 
391   /* Local matrix of all local Schur on subsets transposed */
392   if (!sub_schurs->S_Ej_all) {
393     ierr = MatCreate(PETSC_COMM_SELF,&sub_schurs->S_Ej_all);CHKERRQ(ierr);
394     ierr = MatSetSizes(sub_schurs->S_Ej_all,PETSC_DECIDE,PETSC_DECIDE,local_size,local_size);CHKERRQ(ierr);
395     ierr = MatSetType(sub_schurs->S_Ej_all,MATAIJ);CHKERRQ(ierr);
396     ierr = MatSeqAIJSetPreallocation(sub_schurs->S_Ej_all,0,nnz);CHKERRQ(ierr);
397   } else {
398     ierr = MatZeroEntries(sub_schurs->S_Ej_all);CHKERRQ(ierr);
399   }
400 
401   S_Ej_tilda_all = 0;
402   S_Ej_inv_all = 0;
403   work = NULL;
404   pivots = NULL;
405   if (sub_schurs->n_subs && deluxe) { /* workspace needed only for GETRI */
406     PetscScalar lwork;
407 
408     B_lwork = -1;
409     ierr = PetscBLASIntCast(max_subset_size_seq,&B_N);CHKERRQ(ierr);
410     ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
411     PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,work,&B_N,pivots,&lwork,&B_lwork,&B_ierr));
412     ierr = PetscFPTrapPop();CHKERRQ(ierr);
413     if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to GETRI Lapack routine %d",(int)B_ierr);
414     ierr = PetscBLASIntCast((PetscInt)lwork,&B_lwork);CHKERRQ(ierr);
415     ierr = PetscMalloc2(B_lwork,&work,max_subset_size_seq,&pivots);CHKERRQ(ierr);
416   }
417 
418   ierr = PetscBTMemzero(sub_schurs->n_subs,sub_schurs->computed_Stilda_subs);CHKERRQ(ierr);
419   if (!sub_schurs->use_mumps) {
420     if (sub_schurs->n_subs_seq) {
421       PetscScalar *work;
422       PetscInt    *dummy_idx;
423 
424       /* Work arrays */
425       ierr = PetscMalloc2(max_subset_size_seq,&dummy_idx,max_subset_size_seq*max_subset_size_seq,&work);CHKERRQ(ierr);
426 
427       /* Loop on local problems end compute Schur complements explicitly */
428       local_size = 0;
429       for (i=0;i<sub_schurs->n_subs_seq;i++) {
430         Mat       S_Ej_expl;
431         PetscInt  j,lpi = sub_schurs->index_sequential[i];
432         PetscBool Sdense;
433 
434         ierr = ISGetLocalSize(sub_schurs->is_subs[lpi],&subset_size);CHKERRQ(ierr);
435         ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&S_Ej_expl);CHKERRQ(ierr);
436         /* TODO manage sym flag */
437         ierr = PCBDDCComputeExplicitSchur(sub_schurs->S_Ej[lpi],PETSC_FALSE,MAT_REUSE_MATRIX,&S_Ej_expl);CHKERRQ(ierr);
438         ierr = PetscObjectTypeCompare((PetscObject)S_Ej_expl,MATSEQDENSE,&Sdense);CHKERRQ(ierr);
439         if (Sdense) {
440           for (j=0;j<subset_size;j++) {
441             dummy_idx[j]=local_size+j;
442           }
443           ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
444         } else {
445           SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not yet implemented for sparse matrices");
446         }
447         local_size += subset_size;
448         ierr = MatDestroy(&sub_schurs->S_Ej[lpi]);CHKERRQ(ierr);
449         sub_schurs->S_Ej[lpi] = S_Ej_expl;
450       }
451       ierr = PetscFree2(dummy_idx,work);CHKERRQ(ierr);
452     }
453   } else {
454     PetscInt    *dummy_idx,n_all;
455     PetscScalar *work;
456 
457     if (compute_Stilda) {
458       ierr = MatCreate(PETSC_COMM_SELF,&S_Ej_tilda_all);CHKERRQ(ierr);
459       ierr = MatSetSizes(S_Ej_tilda_all,PETSC_DECIDE,PETSC_DECIDE,local_size,local_size);CHKERRQ(ierr);
460       ierr = MatSetType(S_Ej_tilda_all,MATAIJ);CHKERRQ(ierr);
461       ierr = MatSeqAIJSetPreallocation(S_Ej_tilda_all,0,nnz);CHKERRQ(ierr);
462       if (deluxe) {
463         ierr = MatCreate(PETSC_COMM_SELF,&S_Ej_inv_all);CHKERRQ(ierr);
464         ierr = MatSetSizes(S_Ej_inv_all,PETSC_DECIDE,PETSC_DECIDE,local_size,local_size);CHKERRQ(ierr);
465         ierr = MatSetType(S_Ej_inv_all,MATAIJ);CHKERRQ(ierr);
466         ierr = MatSeqAIJSetPreallocation(S_Ej_inv_all,0,nnz);CHKERRQ(ierr);
467       }
468     }
469 
470     if (S_all) { /* multilevel */
471       ierr = MatGetSize(S_all,&n_all,NULL);CHKERRQ(ierr);
472     } else {
473       n_all = 0;
474     }
475     local_size = 0;
476     /* Work arrays */
477     ierr = PetscMalloc2(max_subset_size_seq,&dummy_idx,max_subset_size_seq*max_subset_size_seq,&work);CHKERRQ(ierr);
478 
479     for (i=0;i<sub_schurs->n_subs;i++) {
480       IS  is_E;
481       PetscScalar *vals;
482       PetscInt j;
483 
484       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
485       for (j=0;j<subset_size;j++) {
486         dummy_idx[j]=local_size+j;
487       }
488       ierr = ISCreateStride(PetscObjectComm((PetscObject)sub_schurs->is_I),subset_size,local_size,1,&is_E);CHKERRQ(ierr);
489       if (sub_schurs->S_Ej[i]) {
490         ierr = MatGetSubMatrix(S_all,is_E,is_E,MAT_REUSE_MATRIX,&sub_schurs->S_Ej[i]);CHKERRQ(ierr);
491       } else {
492         ierr = MatGetSubMatrix(S_all,is_E,is_E,MAT_INITIAL_MATRIX,&sub_schurs->S_Ej[i]);CHKERRQ(ierr);
493       }
494       ierr = MatDenseGetArray(sub_schurs->S_Ej[i],&vals);CHKERRQ(ierr);
495       ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,vals,INSERT_VALUES);CHKERRQ(ierr);
496       ierr = MatDenseRestoreArray(sub_schurs->S_Ej[i],&vals);CHKERRQ(ierr);
497 
498       if (compute_Stilda && ((PetscBTLookup(sub_schurs->is_edge,i) && use_edges) || use_faces)) {
499         ierr = PetscBTSet(sub_schurs->computed_Stilda_subs,i);CHKERRQ(ierr);
500       }
501 
502       if (PetscBTLookup(sub_schurs->computed_Stilda_subs,i)) {
503         Mat Stilda;
504         if (sub_schurs->n_subs == 1) {
505           ierr = PetscObjectReference((PetscObject)sub_schurs->S_Ej[i]);CHKERRQ(ierr);
506           Stilda = sub_schurs->S_Ej[i];
507         } else {
508           KSP ksp;
509           PC  pc;
510           Mat S_EF,S_FE,S_FF,Stilda_impl;
511           IS  is_F;
512           PetscScalar eps=1.e-8;
513           PetscBool chop=PETSC_FALSE;
514 
515           ierr = ISComplement(is_E,0,n_all,&is_F);CHKERRQ(ierr);
516           ierr = MatGetSubMatrix(S_all,is_F,is_F,MAT_INITIAL_MATRIX,&S_FF);CHKERRQ(ierr);
517           ierr = MatGetSubMatrix(S_all,is_F,is_E,MAT_INITIAL_MATRIX,&S_FE);CHKERRQ(ierr);
518           if (sub_schurs->is_hermitian) { /* just need a placeholder for S_EF */
519             ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,n_all-subset_size,work,&S_EF);CHKERRQ(ierr);
520           } else {
521             ierr = MatGetSubMatrix(S_all,is_E,is_F,MAT_INITIAL_MATRIX,&S_EF);CHKERRQ(ierr);
522           }
523           ierr = ISDestroy(&is_F);CHKERRQ(ierr);
524           if (chop) {
525             ierr = MatChop(S_FF,eps);CHKERRQ(ierr);
526             ierr = MatConvert(S_FF,MATAIJ,MAT_REUSE_MATRIX,&S_FF);CHKERRQ(ierr);
527           }
528           ierr = MatCreateSchurComplement(S_FF,S_FF,S_FE,S_EF,sub_schurs->S_Ej[i],&Stilda_impl);CHKERRQ(ierr);
529           ierr = MatDestroy(&S_FF);CHKERRQ(ierr);
530           ierr = MatDestroy(&S_FE);CHKERRQ(ierr);
531           ierr = MatDestroy(&S_EF);CHKERRQ(ierr);
532           ierr = MatSchurComplementGetKSP(Stilda_impl,&ksp);CHKERRQ(ierr);
533           ierr = KSPSetType(ksp,KSPPREONLY);CHKERRQ(ierr);
534           ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr);
535           if (sub_schurs->is_hermitian && sub_schurs->is_posdef) {
536             ierr = PCSetType(pc,PCCHOLESKY);CHKERRQ(ierr);
537           } else {
538             ierr = PCSetType(pc,PCLU);CHKERRQ(ierr);
539           }
540           if (!chop) {
541             ierr = PCFactorSetUseInPlace(pc,PETSC_TRUE);CHKERRQ(ierr);
542           } else {
543             ierr = PCFactorSetMatSolverPackage(pc,MATSOLVERMUMPS);CHKERRQ(ierr);
544           }
545           ierr = KSPSetUp(ksp);CHKERRQ(ierr);
546           ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&Stilda);CHKERRQ(ierr);
547           ierr = PCBDDCComputeExplicitSchur(Stilda_impl,sub_schurs->is_hermitian,MAT_REUSE_MATRIX,&Stilda);CHKERRQ(ierr);
548           ierr = MatDestroy(&Stilda_impl);CHKERRQ(ierr);
549         }
550 /*
551         PetscViewerPushFormat(PETSC_VIEWER_STDOUT_SELF,PETSC_VIEWER_ASCII_MATLAB);
552         ierr = MatView(Stilda,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
553         PetscViewerPopFormat(PETSC_VIEWER_STDOUT_SELF);
554 */
555         ierr = MatDenseGetArray(Stilda,&vals);CHKERRQ(ierr);
556         if (deluxe) { /* when using deluxe scaling, we need (S_1^-1+S_2^-1)^-1 */
557           PetscScalar *vals2;
558 
559           ierr = MatDenseGetArray(sub_schurs->S_Ej[i],&vals2);CHKERRQ(ierr);
560           /* need to be optimized (cholesky) */
561           ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr);
562           ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
563           if (invert_Stildas) { /* Stildas can be singular */
564             PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,vals,&B_N,pivots,&B_ierr));
565             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
566             PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,vals,&B_N,pivots,work,&B_lwork,&B_ierr));
567             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
568           }
569           PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,vals2,&B_N,pivots,&B_ierr));
570           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
571           PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,vals2,&B_N,pivots,work,&B_lwork,&B_ierr));
572           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
573           ierr = PetscFPTrapPop();CHKERRQ(ierr);
574           ierr = MatSetValues(S_Ej_inv_all,subset_size,dummy_idx,subset_size,dummy_idx,vals2,INSERT_VALUES);CHKERRQ(ierr);
575           ierr = MatDenseRestoreArray(sub_schurs->S_Ej[i],&vals2);CHKERRQ(ierr);
576         }
577         ierr = MatSetValues(S_Ej_tilda_all,subset_size,dummy_idx,subset_size,dummy_idx,vals,INSERT_VALUES);CHKERRQ(ierr);
578         ierr = MatDenseRestoreArray(Stilda,&vals);CHKERRQ(ierr);
579         ierr = MatDestroy(&Stilda);CHKERRQ(ierr);
580       }
581       ierr = ISDestroy(&is_E);CHKERRQ(ierr);
582       local_size += subset_size;
583     }
584     ierr = PetscFree2(dummy_idx,work);CHKERRQ(ierr);
585   }
586   ierr = PetscFree(nnz);CHKERRQ(ierr);
587   ierr = MatDestroy(&S_all);CHKERRQ(ierr);
588   ierr = MatAssemblyBegin(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
589   ierr = MatAssemblyEnd(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
590   if (S_Ej_tilda_all) {
591     ierr = MatAssemblyBegin(S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
592     ierr = MatAssemblyEnd(S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
593   }
594   if (S_Ej_inv_all) {
595     ierr = MatAssemblyBegin(S_Ej_inv_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
596     ierr = MatAssemblyEnd(S_Ej_inv_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
597   }
598 
599   /* Global matrix of all assembled Schur on subsets */
600   ierr = PCBDDCSubsetNumbering(PetscObjectComm((PetscObject)sub_schurs->l2gmap),sub_schurs->l2gmap,local_size,all_local_idx_N+extra,PETSC_NULL,&global_size,&all_local_idx_G);CHKERRQ(ierr);
601   ierr = ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)sub_schurs->l2gmap),1,local_size,all_local_idx_G,PETSC_COPY_VALUES,&l2gmap_subsets);CHKERRQ(ierr);
602   ierr = MatCreateIS(PetscObjectComm((PetscObject)sub_schurs->l2gmap),1,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size,l2gmap_subsets,&work_mat);CHKERRQ(ierr);
603   ierr = ISLocalToGlobalMappingDestroy(&l2gmap_subsets);CHKERRQ(ierr);
604   ierr = MatISSetLocalMat(work_mat,sub_schurs->S_Ej_all);CHKERRQ(ierr);
605   ierr = MatCreate(PetscObjectComm((PetscObject)work_mat),&global_schur_subsets);CHKERRQ(ierr);
606   ierr = MatSetSizes(global_schur_subsets,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size);CHKERRQ(ierr);
607   ierr = MatSetType(global_schur_subsets,MATMPIAIJ);CHKERRQ(ierr);
608   ierr = MatISSetMPIXAIJPreallocation_Private(work_mat,global_schur_subsets,PETSC_TRUE);CHKERRQ(ierr);
609   ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
610   /* Get local part of (\sum_j S_Ej) */
611   ierr = PetscFree(all_local_idx_N);CHKERRQ(ierr);
612   ierr = ISCreateGeneral(PetscObjectComm((PetscObject)sub_schurs->l2gmap),local_size,all_local_idx_G,PETSC_OWN_POINTER,&temp_is);CHKERRQ(ierr);
613   if (sub_schurs->sum_S_Ej_all) {
614     ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,MAT_REUSE_MATRIX,&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
615   } else {
616     ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,MAT_INITIAL_MATRIX,&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
617   }
618 
619   if (S_Ej_tilda_all) {
620     ierr = MatISSetLocalMat(work_mat,S_Ej_tilda_all);CHKERRQ(ierr);
621     ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
622     if (sub_schurs->sum_S_Ej_tilda_all) {
623       ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,MAT_REUSE_MATRIX,&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
624     } else {
625       ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,MAT_INITIAL_MATRIX,&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
626     }
627   }
628   if (S_Ej_inv_all) {
629     PetscInt    cum;
630     PetscScalar *array,*array2;
631     ierr = MatISSetLocalMat(work_mat,S_Ej_inv_all);CHKERRQ(ierr);
632     ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
633     if (sub_schurs->sum_S_Ej_inv_all) {
634       ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,MAT_REUSE_MATRIX,&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
635     } else {
636       ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,MAT_INITIAL_MATRIX,&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
637     }
638     /* invert blocks of sum_S_Ej_inv_all */
639     ierr = MatSeqAIJGetArray(sub_schurs->sum_S_Ej_inv_all,&array);CHKERRQ(ierr);
640     ierr = MatSeqAIJGetArray(sub_schurs->sum_S_Ej_tilda_all,&array2);CHKERRQ(ierr);
641     cum = 0;
642     for (i=0;i<sub_schurs->n_subs;i++) {
643       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
644       if (PetscBTLookup(sub_schurs->computed_Stilda_subs,i)) {
645         /* need to be optimized (cholesky) */
646         ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr);
647         ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
648         PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,array+cum,&B_N,pivots,&B_ierr));
649         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
650         PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,array+cum,&B_N,pivots,work,&B_lwork,&B_ierr));
651         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
652         if (invert_Stildas) { /* Stildas can be singular */
653           PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,array2+cum,&B_N,pivots,&B_ierr));
654           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
655           PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,array2+cum,&B_N,pivots,work,&B_lwork,&B_ierr));
656           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
657         }
658         ierr = PetscFPTrapPop();CHKERRQ(ierr);
659       }
660       cum += subset_size*subset_size;
661     }
662     ierr = MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_inv_all,&array);CHKERRQ(ierr);
663     ierr = MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_tilda_all,&array2);CHKERRQ(ierr);
664   }
665 
666   ierr = PetscFree2(work,pivots);CHKERRQ(ierr);
667   ierr = MatDestroy(&global_schur_subsets);CHKERRQ(ierr);
668   ierr = MatDestroy(&S_Ej_tilda_all);CHKERRQ(ierr);
669   ierr = MatDestroy(&S_Ej_inv_all);CHKERRQ(ierr);
670   ierr = MatDestroy(&work_mat);CHKERRQ(ierr);
671   ierr = ISDestroy(&temp_is);CHKERRQ(ierr);
672   PetscFunctionReturn(0);
673 }
674 
675 #undef __FUNCT__
676 #define __FUNCT__ "PCBDDCSubSchursInit"
677 PetscErrorCode PCBDDCSubSchursInit(PCBDDCSubSchurs sub_schurs, Mat A, Mat S, IS is_I, IS is_B, PCBDDCGraph graph, ISLocalToGlobalMapping BtoNmap, PetscInt seqthreshold)
678 {
679   IS                  *faces,*edges,*all_cc,vertices;
680   PetscInt            *index_sequential,*index_parallel;
681   PetscInt            *auxlocal_sequential,*auxlocal_parallel;
682   PetscInt            *auxglobal_sequential,*auxglobal_parallel;
683   PetscInt            *auxmapping;
684   PetscInt            i,max_subset_size;
685   PetscInt            n_sequential_problems,n_local_sequential_problems,n_parallel_problems,n_local_parallel_problems;
686   PetscInt            n_faces,n_edges,n_all_cc;
687   PetscBool           is_sorted;
688   PetscErrorCode  ierr;
689 
690   PetscFunctionBegin;
691   ierr = ISSorted(is_I,&is_sorted);CHKERRQ(ierr);
692   if (!is_sorted) {
693     SETERRQ(PetscObjectComm((PetscObject)is_I),PETSC_ERR_PLIB,"IS for I dofs should be shorted");
694   }
695   ierr = ISSorted(is_B,&is_sorted);CHKERRQ(ierr);
696   if (!is_sorted) {
697     SETERRQ(PetscObjectComm((PetscObject)is_B),PETSC_ERR_PLIB,"IS for B dofs should be shorted");
698   }
699 
700   /* reset any previous data */
701   ierr = PCBDDCSubSchursReset(sub_schurs);CHKERRQ(ierr);
702 
703   /* get index sets for faces and edges */
704   ierr = PCBDDCGraphGetCandidatesIS(graph,&n_faces,&faces,&n_edges,&edges,&vertices);CHKERRQ(ierr);
705   n_all_cc = n_faces+n_edges;
706   ierr = PetscBTCreate(n_all_cc,&sub_schurs->is_edge);CHKERRQ(ierr);
707   ierr = PetscBTCreate(n_all_cc,&sub_schurs->computed_Stilda_subs);CHKERRQ(ierr);
708   ierr = PetscMalloc1(n_all_cc,&all_cc);CHKERRQ(ierr);
709   for (i=0;i<n_faces;i++) {
710     all_cc[i] = faces[i];
711   }
712   for (i=0;i<n_edges;i++) {
713     all_cc[n_faces+i] = edges[i];
714     ierr = PetscBTSet(sub_schurs->is_edge,n_faces+i);CHKERRQ(ierr);
715   }
716   ierr = PetscFree(faces);CHKERRQ(ierr);
717   ierr = PetscFree(edges);CHKERRQ(ierr);
718 
719   /* map interface's subsets */
720   max_subset_size = 0;
721   for (i=0;i<n_all_cc;i++) {
722     PetscInt subset_size;
723     ierr = ISGetLocalSize(all_cc[i],&subset_size);CHKERRQ(ierr);
724     max_subset_size = PetscMax(max_subset_size,subset_size);
725   }
726   ierr = PetscMalloc1(max_subset_size,&auxmapping);CHKERRQ(ierr);
727   ierr = PetscMalloc2(graph->ncc,&auxlocal_sequential,
728                       graph->ncc,&auxlocal_parallel);CHKERRQ(ierr);
729   ierr = PetscMalloc2(graph->ncc,&index_sequential,
730                       graph->ncc,&index_parallel);CHKERRQ(ierr);
731 
732   /* if threshold is negative uses all sequential problems (possibly using MUMPS) */
733   sub_schurs->use_mumps = PETSC_FALSE;
734   if (seqthreshold < 0) {
735     seqthreshold = max_subset_size;
736 #if defined(PETSC_HAVE_MUMPS)
737     sub_schurs->use_mumps = !!A;
738 #endif
739   }
740 
741 
742   /* determine which problem has to be solved in parallel or sequentially */
743   n_local_sequential_problems = 0;
744   n_local_parallel_problems = 0;
745   for (i=0;i<n_all_cc;i++) {
746     PetscInt       subset_size,j,min_loc = 0;
747     const PetscInt *idxs;
748 
749     ierr = ISGetLocalSize(all_cc[i],&subset_size);CHKERRQ(ierr);
750     ierr = ISGetIndices(all_cc[i],&idxs);CHKERRQ(ierr);
751     ierr = ISLocalToGlobalMappingApply(graph->l2gmap,subset_size,idxs,auxmapping);CHKERRQ(ierr);
752     for (j=1;j<subset_size;j++) {
753       if (auxmapping[j]<auxmapping[min_loc]) {
754         min_loc = j;
755       }
756     }
757     if (subset_size > seqthreshold) {
758       index_parallel[n_local_parallel_problems] = i;
759       auxlocal_parallel[n_local_parallel_problems] = idxs[min_loc];
760       n_local_parallel_problems++;
761     } else {
762       index_sequential[n_local_sequential_problems] = i;
763       auxlocal_sequential[n_local_sequential_problems] = idxs[min_loc];
764       n_local_sequential_problems++;
765     }
766     ierr = ISRestoreIndices(all_cc[i],&idxs);CHKERRQ(ierr);
767   }
768 
769   /* Number parallel problems */
770   auxglobal_parallel = 0;
771   ierr = PCBDDCSubsetNumbering(PetscObjectComm((PetscObject)graph->l2gmap),graph->l2gmap,n_local_parallel_problems,auxlocal_parallel,PETSC_NULL,&n_parallel_problems,&auxglobal_parallel);CHKERRQ(ierr);
772 
773   /* Number sequential problems */
774   auxglobal_sequential = 0;
775   ierr = PCBDDCSubsetNumbering(PetscObjectComm((PetscObject)graph->l2gmap),graph->l2gmap,n_local_sequential_problems,auxlocal_sequential,PETSC_NULL,&n_sequential_problems,&auxglobal_sequential);CHKERRQ(ierr);
776 
777   /* update info in sub_schurs */
778   if (sub_schurs->use_mumps && A) {
779     PetscBool isseqaij;
780 
781     ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
782     if (isseqaij) {
783       ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
784       sub_schurs->A = A;
785     } else { /* SeqBAIJ matrices does not support symmetry checking, SeqSBAIJ does not support MatPermute */
786       ierr = MatConvert(A,MATSEQAIJ,MAT_INITIAL_MATRIX,&sub_schurs->A);CHKERRQ(ierr);
787     }
788   }
789   ierr = PetscObjectReference((PetscObject)S);CHKERRQ(ierr);
790   sub_schurs->S = S;
791   ierr = PetscObjectReference((PetscObject)is_I);CHKERRQ(ierr);
792   sub_schurs->is_I = is_I;
793   ierr = PetscObjectReference((PetscObject)is_B);CHKERRQ(ierr);
794   sub_schurs->is_B = is_B;
795   ierr = PetscObjectReference((PetscObject)graph->l2gmap);CHKERRQ(ierr);
796   sub_schurs->l2gmap = graph->l2gmap;
797   ierr = PetscObjectReference((PetscObject)BtoNmap);CHKERRQ(ierr);
798   sub_schurs->BtoNmap = BtoNmap;
799   sub_schurs->n_subs_seq = n_local_sequential_problems;
800   sub_schurs->n_subs_par = n_local_parallel_problems;
801   sub_schurs->n_subs_seq_g = n_sequential_problems;
802   sub_schurs->n_subs_par_g = n_parallel_problems;
803   sub_schurs->n_subs = sub_schurs->n_subs_seq + sub_schurs->n_subs_par;
804   sub_schurs->is_subs = all_cc;
805   if (!sub_schurs->use_mumps) { /* for adaptive selection */
806     for (i=0;i<sub_schurs->n_subs;i++) {
807       ierr = ISSort(sub_schurs->is_subs[i]);CHKERRQ(ierr);
808     }
809   }
810   sub_schurs->is_Ej_com = vertices;
811   sub_schurs->index_sequential = index_sequential;
812   sub_schurs->index_parallel = index_parallel;
813   sub_schurs->auxglobal_sequential = auxglobal_sequential;
814   sub_schurs->auxglobal_parallel = auxglobal_parallel;
815 
816 
817   /* allocate space for schur complements */
818   ierr = PetscCalloc1(sub_schurs->n_subs,&sub_schurs->S_Ej);CHKERRQ(ierr);
819   sub_schurs->S_Ej_all = NULL;
820   sub_schurs->sum_S_Ej_all = NULL;
821   sub_schurs->sum_S_Ej_inv_all = NULL;
822   sub_schurs->sum_S_Ej_tilda_all = NULL;
823   sub_schurs->is_Ej_all = NULL;
824 
825   /* free workspace */
826   ierr = PetscFree(auxmapping);CHKERRQ(ierr);
827   ierr = PetscFree2(auxlocal_sequential,auxlocal_parallel);CHKERRQ(ierr);
828 
829   PetscFunctionReturn(0);
830 }
831 
832 #undef __FUNCT__
833 #define __FUNCT__ "PCBDDCSubSchursCreate"
834 PetscErrorCode PCBDDCSubSchursCreate(PCBDDCSubSchurs *sub_schurs)
835 {
836   PCBDDCSubSchurs schurs_ctx;
837   PetscErrorCode  ierr;
838 
839   PetscFunctionBegin;
840   ierr = PetscNew(&schurs_ctx);CHKERRQ(ierr);
841   schurs_ctx->n_subs = 0;
842   *sub_schurs = schurs_ctx;
843   PetscFunctionReturn(0);
844 }
845 
846 #undef __FUNCT__
847 #define __FUNCT__ "PCBDDCSubSchursDestroy"
848 PetscErrorCode PCBDDCSubSchursDestroy(PCBDDCSubSchurs *sub_schurs)
849 {
850   PetscErrorCode ierr;
851 
852   PetscFunctionBegin;
853   ierr = PCBDDCSubSchursReset(*sub_schurs);CHKERRQ(ierr);
854   ierr = PetscFree(*sub_schurs);CHKERRQ(ierr);
855   PetscFunctionReturn(0);
856 }
857 
858 #undef __FUNCT__
859 #define __FUNCT__ "PCBDDCSubSchursReset"
860 PetscErrorCode PCBDDCSubSchursReset(PCBDDCSubSchurs sub_schurs)
861 {
862   PetscInt       i;
863   PetscErrorCode ierr;
864 
865   PetscFunctionBegin;
866   ierr = MatDestroy(&sub_schurs->A);CHKERRQ(ierr);
867   ierr = MatDestroy(&sub_schurs->S);CHKERRQ(ierr);
868   ierr = ISDestroy(&sub_schurs->is_I);CHKERRQ(ierr);
869   ierr = ISDestroy(&sub_schurs->is_B);CHKERRQ(ierr);
870   ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->l2gmap);CHKERRQ(ierr);
871   ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->BtoNmap);CHKERRQ(ierr);
872   ierr = MatDestroy(&sub_schurs->S_Ej_all);CHKERRQ(ierr);
873   ierr = MatDestroy(&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
874   ierr = MatDestroy(&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
875   ierr = MatDestroy(&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
876   ierr = ISDestroy(&sub_schurs->is_Ej_all);CHKERRQ(ierr);
877   ierr = ISDestroy(&sub_schurs->is_Ej_com);CHKERRQ(ierr);
878   ierr = PetscBTDestroy(&sub_schurs->is_edge);CHKERRQ(ierr);
879   ierr = PetscBTDestroy(&sub_schurs->computed_Stilda_subs);CHKERRQ(ierr);
880   for (i=0;i<sub_schurs->n_subs;i++) {
881     ierr = ISDestroy(&sub_schurs->is_subs[i]);CHKERRQ(ierr);
882     ierr = MatDestroy(&sub_schurs->S_Ej[i]);CHKERRQ(ierr);
883   }
884   ierr = ISDestroy(&sub_schurs->is_I_layer);CHKERRQ(ierr);
885   if (sub_schurs->n_subs) {
886     ierr = PetscFree(sub_schurs->is_subs);CHKERRQ(ierr);
887     ierr = PetscFree(sub_schurs->S_Ej);CHKERRQ(ierr);
888   }
889   ierr = PetscFree2(sub_schurs->index_sequential,sub_schurs->index_parallel);CHKERRQ(ierr);
890   ierr = PetscFree(sub_schurs->auxglobal_sequential);CHKERRQ(ierr);
891   ierr = PetscFree(sub_schurs->auxglobal_parallel);CHKERRQ(ierr);
892   sub_schurs->n_subs = 0;
893   PetscFunctionReturn(0);
894 }
895 
896 #undef __FUNCT__
897 #define __FUNCT__ "PCBDDCAdjGetNextLayer_Private"
898 PETSC_STATIC_INLINE PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt* queue_tip,PetscInt n_prev,PetscBT touched,PetscInt* xadj,PetscInt* adjncy,PetscInt* n_added)
899 {
900   PetscInt       i,j,n;
901   PetscErrorCode ierr;
902 
903   PetscFunctionBegin;
904   n = 0;
905   for (i=-n_prev;i<0;i++) {
906     PetscInt start_dof = queue_tip[i];
907     for (j=xadj[start_dof];j<xadj[start_dof+1];j++) {
908       PetscInt dof = adjncy[j];
909       if (!PetscBTLookup(touched,dof)) {
910         ierr = PetscBTSet(touched,dof);CHKERRQ(ierr);
911         queue_tip[n] = dof;
912         n++;
913       }
914     }
915   }
916   *n_added = n;
917   PetscFunctionReturn(0);
918 }
919