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