xref: /petsc/src/ksp/pc/impls/bddc/bddcschurs.c (revision 06a4b1fae8f1b3303b6dfb6f4a5d2fdf4f3d2719)
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            *Bwork;
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 
276       /* IE block */
277       for (i=0;i<sub_schurs->n_subs;i++) {
278         ierr = MatGetSubMatrix(A_IB,is_I,is_subset_B[i],MAT_INITIAL_MATRIX,&AE_IE[i]);CHKERRQ(ierr);
279       }
280 
281       /* EI block */
282       for (i=0;i<sub_schurs->n_subs;i++) {
283         if (sub_schurs->is_hermitian) {
284           ierr = MatCreateTranspose(AE_IE[i],&AE_EI[i]);CHKERRQ(ierr);
285         } else {
286           ierr = MatGetSubMatrix(A_BI,is_subset_B[i],is_I,MAT_INITIAL_MATRIX,&AE_EI[i]);CHKERRQ(ierr);
287         }
288       }
289 
290       /* setup Schur complements on subset */
291       for (i=0;i<sub_schurs->n_subs;i++) {
292         ierr = MatDestroy(&sub_schurs->S_Ej[i]);CHKERRQ(ierr);
293         ierr = MatCreateSchurComplement(AE_II,AE_II,AE_IE[i],AE_EI[i],AE_EE[i],&sub_schurs->S_Ej[i]);CHKERRQ(ierr);
294         ierr = MatDestroy(&AE_EE[i]);CHKERRQ(ierr);
295         ierr = MatDestroy(&AE_IE[i]);CHKERRQ(ierr);
296         ierr = MatDestroy(&AE_EI[i]);CHKERRQ(ierr);
297         if (AE_II == A_II) { /* we can reuse the same ksp */
298           KSP ksp;
299           ierr = MatSchurComplementGetKSP(sub_schurs->S,&ksp);CHKERRQ(ierr);
300           ierr = MatSchurComplementSetKSP(sub_schurs->S_Ej[i],ksp);CHKERRQ(ierr);
301         } else { /* build new ksp object which inherits ksp and pc types from the original one */
302           KSP      origksp,schurksp;
303           PC       origpc,schurpc;
304           KSPType  ksp_type;
305           PCType   pc_type;
306           PetscInt n_internal;
307 
308           ierr = MatSchurComplementGetKSP(sub_schurs->S,&origksp);CHKERRQ(ierr);
309           ierr = MatSchurComplementGetKSP(sub_schurs->S_Ej[i],&schurksp);CHKERRQ(ierr);
310           ierr = KSPGetType(origksp,&ksp_type);CHKERRQ(ierr);
311           ierr = KSPSetType(schurksp,ksp_type);CHKERRQ(ierr);
312           ierr = KSPGetPC(schurksp,&schurpc);CHKERRQ(ierr);
313           ierr = KSPGetPC(origksp,&origpc);CHKERRQ(ierr);
314           ierr = PCGetType(origpc,&pc_type);CHKERRQ(ierr);
315           ierr = PCSetType(schurpc,pc_type);CHKERRQ(ierr);
316           ierr = ISGetSize(is_I,&n_internal);CHKERRQ(ierr);
317           if (n_internal) { /* UMFPACK gives error with 0 sized problems */
318             MatSolverPackage solver=NULL;
319             ierr = PCFactorGetMatSolverPackage(origpc,(const MatSolverPackage*)&solver);CHKERRQ(ierr);
320             if (solver) {
321               ierr = PCFactorSetMatSolverPackage(schurpc,solver);CHKERRQ(ierr);
322             }
323           }
324           ierr = KSPSetUp(schurksp);CHKERRQ(ierr);
325         }
326       }
327       /* free */
328       ierr = ISDestroy(&is_I);CHKERRQ(ierr);
329       ierr = MatDestroy(&AE_II);CHKERRQ(ierr);
330       for (i=0;i<sub_schurs->n_subs;i++) {
331         ierr = ISDestroy(&is_subset_B[i]);CHKERRQ(ierr);
332       }
333       ierr = PetscFree4(is_subset_B,AE_IE,AE_EI,AE_EE);CHKERRQ(ierr);
334 #if defined(PETSC_HAVE_MUMPS)
335     } else {
336       Mat           A,F;
337       IS            is_A_all;
338       PetscInt      *idxs_schur,n_I;
339 
340       /* get working mat */
341       ierr = ISGetLocalSize(sub_schurs->is_I_layer,&n_I);CHKERRQ(ierr);
342       ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size+n_I,all_local_idx_N,PETSC_COPY_VALUES,&is_A_all);CHKERRQ(ierr);
343       ierr = MatGetSubMatrixUnsorted(sub_schurs->A,is_A_all,is_A_all,MAT_INITIAL_MATRIX,&A);CHKERRQ(ierr);
344       ierr = ISDestroy(&is_A_all);CHKERRQ(ierr);
345 
346       if (n_I) {
347         if (sub_schurs->is_hermitian && sub_schurs->is_posdef) {
348           ierr = MatGetFactor(A,MATSOLVERMUMPS,MAT_FACTOR_CHOLESKY,&F);CHKERRQ(ierr);
349         } else {
350           ierr = MatGetFactor(A,MATSOLVERMUMPS,MAT_FACTOR_LU,&F);CHKERRQ(ierr);
351         }
352 
353         /* subsets ordered last */
354         ierr = PetscMalloc1(local_size,&idxs_schur);CHKERRQ(ierr);
355         for (i=0;i<local_size;i++) {
356           idxs_schur[i] = n_I+i+1;
357         }
358         ierr = MatMumpsSetSchurIndices(F,local_size,idxs_schur);CHKERRQ(ierr);
359         ierr = PetscFree(idxs_schur);CHKERRQ(ierr);
360 
361         /* factorization step */
362         if (sub_schurs->is_hermitian && sub_schurs->is_posdef) {
363           ierr = MatCholeskyFactorSymbolic(F,A,NULL,NULL);CHKERRQ(ierr);
364           ierr = MatCholeskyFactorNumeric(F,A,NULL);CHKERRQ(ierr);
365         } else {
366           ierr = MatLUFactorSymbolic(F,A,NULL,NULL,NULL);CHKERRQ(ierr);
367           ierr = MatLUFactorNumeric(F,A,NULL);CHKERRQ(ierr);
368         }
369 
370         /* get explicit Schur Complement computed during numeric factorization */
371         ierr = MatMumpsGetSchurComplement(F,&S_all);CHKERRQ(ierr);
372         ierr = MatDestroy(&F);CHKERRQ(ierr);
373       } else {
374         ierr = MatConvert(A,MATSEQDENSE,MAT_INITIAL_MATRIX,&S_all);CHKERRQ(ierr);
375       }
376       ierr = MatDestroy(&A);CHKERRQ(ierr);
377 #endif
378     }
379   } else {
380     ierr = PetscFree(nnz);CHKERRQ(ierr);
381     ierr = PetscFree(all_local_idx_N);CHKERRQ(ierr);
382   }
383 
384   if (!sub_schurs->n_subs_seq_g) {
385     PetscFunctionReturn(0);
386   }
387 
388   /* subset indices in local boundary numbering */
389   if (!sub_schurs->is_Ej_all) {
390     ierr = ISGlobalToLocalMappingApply(sub_schurs->BtoNmap,IS_GTOLM_DROP,local_size,all_local_idx_N+extra,&subset_size,all_local_idx_B);CHKERRQ(ierr);
391     if (subset_size != local_size) {
392       SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in sub_schurs serial (BtoNmap)! %d != %d\n",subset_size,local_size);
393     }
394     ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size,all_local_idx_B,PETSC_OWN_POINTER,&sub_schurs->is_Ej_all);CHKERRQ(ierr);
395   }
396 
397   /* Local matrix of all local Schur on subsets transposed */
398   if (!sub_schurs->S_Ej_all) {
399     ierr = MatCreate(PETSC_COMM_SELF,&sub_schurs->S_Ej_all);CHKERRQ(ierr);
400     ierr = MatSetSizes(sub_schurs->S_Ej_all,PETSC_DECIDE,PETSC_DECIDE,local_size,local_size);CHKERRQ(ierr);
401     ierr = MatSetType(sub_schurs->S_Ej_all,MATAIJ);CHKERRQ(ierr);
402     ierr = MatSeqAIJSetPreallocation(sub_schurs->S_Ej_all,0,nnz);CHKERRQ(ierr);
403   } else {
404     ierr = MatZeroEntries(sub_schurs->S_Ej_all);CHKERRQ(ierr);
405   }
406 
407   S_Ej_tilda_all = 0;
408   S_Ej_inv_all = 0;
409   Bwork = NULL;
410   pivots = NULL;
411   if (sub_schurs->n_subs && deluxe && compute_Stilda) { /* workspace needed only for GETRI */
412     PetscScalar lwork;
413 
414     B_lwork = -1;
415     ierr = PetscBLASIntCast(max_subset_size_seq,&B_N);CHKERRQ(ierr);
416     ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
417     PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,Bwork,&B_N,pivots,&lwork,&B_lwork,&B_ierr));
418     ierr = PetscFPTrapPop();CHKERRQ(ierr);
419     if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to GETRI Lapack routine %d",(int)B_ierr);
420     ierr = PetscBLASIntCast((PetscInt)lwork,&B_lwork);CHKERRQ(ierr);
421     ierr = PetscMalloc2(B_lwork,&Bwork,max_subset_size_seq,&pivots);CHKERRQ(ierr);
422   }
423 
424   ierr = PetscBTMemzero(sub_schurs->n_subs,sub_schurs->computed_Stilda_subs);CHKERRQ(ierr);
425   if (!sub_schurs->use_mumps) {
426     if (sub_schurs->n_subs_seq) {
427       PetscScalar *work;
428       PetscInt    *dummy_idx;
429 
430       /* Work arrays */
431       ierr = PetscMalloc2(max_subset_size_seq,&dummy_idx,max_subset_size_seq*max_subset_size_seq,&work);CHKERRQ(ierr);
432 
433       /* Loop on local problems end compute Schur complements explicitly */
434       local_size = 0;
435       for (i=0;i<sub_schurs->n_subs_seq;i++) {
436         Mat       S_Ej_expl;
437         PetscInt  j,lpi = sub_schurs->index_sequential[i];
438         PetscBool Sdense;
439 
440         ierr = ISGetLocalSize(sub_schurs->is_subs[lpi],&subset_size);CHKERRQ(ierr);
441         ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&S_Ej_expl);CHKERRQ(ierr);
442         ierr = PCBDDCComputeExplicitSchur(sub_schurs->S_Ej[lpi],sub_schurs->is_hermitian,MAT_REUSE_MATRIX,&S_Ej_expl);CHKERRQ(ierr);
443         ierr = PetscObjectTypeCompare((PetscObject)S_Ej_expl,MATSEQDENSE,&Sdense);CHKERRQ(ierr);
444         if (Sdense) {
445           for (j=0;j<subset_size;j++) {
446             dummy_idx[j]=local_size+j;
447           }
448           ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
449         } else {
450           SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not yet implemented for sparse matrices");
451         }
452         local_size += subset_size;
453         ierr = MatDestroy(&sub_schurs->S_Ej[lpi]);CHKERRQ(ierr);
454         sub_schurs->S_Ej[lpi] = S_Ej_expl;
455       }
456       ierr = PetscFree2(dummy_idx,work);CHKERRQ(ierr);
457     }
458   } else {
459     PetscInt    *dummy_idx,n_all;
460     PetscScalar *work;
461 
462     if (compute_Stilda) {
463       ierr = MatCreate(PETSC_COMM_SELF,&S_Ej_tilda_all);CHKERRQ(ierr);
464       ierr = MatSetSizes(S_Ej_tilda_all,PETSC_DECIDE,PETSC_DECIDE,local_size,local_size);CHKERRQ(ierr);
465       ierr = MatSetType(S_Ej_tilda_all,MATAIJ);CHKERRQ(ierr);
466       ierr = MatSeqAIJSetPreallocation(S_Ej_tilda_all,0,nnz);CHKERRQ(ierr);
467       if (deluxe) {
468         ierr = MatCreate(PETSC_COMM_SELF,&S_Ej_inv_all);CHKERRQ(ierr);
469         ierr = MatSetSizes(S_Ej_inv_all,PETSC_DECIDE,PETSC_DECIDE,local_size,local_size);CHKERRQ(ierr);
470         ierr = MatSetType(S_Ej_inv_all,MATAIJ);CHKERRQ(ierr);
471         ierr = MatSeqAIJSetPreallocation(S_Ej_inv_all,0,nnz);CHKERRQ(ierr);
472       }
473     }
474 
475     if (S_all) { /* multilevel */
476       ierr = MatGetSize(S_all,&n_all,NULL);CHKERRQ(ierr);
477     } else {
478       n_all = 0;
479     }
480     local_size = 0;
481 
482     /* Work arrays */
483     if (compute_Stilda) {
484       ierr = PetscMalloc2(max_subset_size_seq,&dummy_idx,n_all*n_all,&work);CHKERRQ(ierr);
485     } else {
486       ierr = PetscMalloc2(max_subset_size_seq,&dummy_idx,0,&work);CHKERRQ(ierr);
487     }
488 
489     for (i=0;i<sub_schurs->n_subs;i++) {
490       IS  is_E;
491       PetscScalar *vals;
492       PetscInt j;
493 
494       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
495       for (j=0;j<subset_size;j++) {
496         dummy_idx[j]=local_size+j;
497       }
498       ierr = ISCreateStride(PetscObjectComm((PetscObject)sub_schurs->is_I),subset_size,local_size,1,&is_E);CHKERRQ(ierr);
499       if (sub_schurs->S_Ej[i]) {
500         ierr = MatGetSubMatrix(S_all,is_E,is_E,MAT_REUSE_MATRIX,&sub_schurs->S_Ej[i]);CHKERRQ(ierr);
501       } else {
502         ierr = MatGetSubMatrix(S_all,is_E,is_E,MAT_INITIAL_MATRIX,&sub_schurs->S_Ej[i]);CHKERRQ(ierr);
503       }
504       ierr = MatDenseGetArray(sub_schurs->S_Ej[i],&vals);CHKERRQ(ierr);
505       ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,vals,INSERT_VALUES);CHKERRQ(ierr);
506       ierr = MatDenseRestoreArray(sub_schurs->S_Ej[i],&vals);CHKERRQ(ierr);
507 
508       if (compute_Stilda && ((PetscBTLookup(sub_schurs->is_edge,i) && use_edges) || use_faces)) {
509         Mat Stilda;
510 
511         ierr = PetscBTSet(sub_schurs->computed_Stilda_subs,i);CHKERRQ(ierr);
512         if (sub_schurs->n_subs == 1) {
513           ierr = PetscObjectReference((PetscObject)sub_schurs->S_Ej[i]);CHKERRQ(ierr);
514           Stilda = sub_schurs->S_Ej[i];
515         } else {
516           KSP         ksp;
517           PC          pc;
518           Mat         S_EF,S_FE,S_FF,Stilda_impl;
519           IS          is_F;
520           PetscScalar eps=1.e-8;
521           PetscInt    nc,cum;
522           PetscBool   chop=PETSC_FALSE;
523 
524           ierr = ISComplement(is_E,0,n_all,&is_F);CHKERRQ(ierr);
525           nc   = n_all - subset_size;
526           cum  = 0;
527           ierr = MatCreateSeqDense(PETSC_COMM_SELF,nc,nc,work+cum,&S_FF);CHKERRQ(ierr);
528           cum  = nc*nc;
529           ierr = MatCreateSeqDense(PETSC_COMM_SELF,nc,subset_size,work+cum,&S_FE);CHKERRQ(ierr);
530           cum  = nc*subset_size;
531           ierr = MatGetSubMatrix(S_all,is_F,is_F,MAT_REUSE_MATRIX,&S_FF);CHKERRQ(ierr);
532           ierr = MatGetSubMatrix(S_all,is_F,is_E,MAT_REUSE_MATRIX,&S_FE);CHKERRQ(ierr);
533           if (sub_schurs->is_hermitian) { /* just need a placeholder for S_EF */
534             ierr = MatCreateTranspose(S_FE,&S_EF);CHKERRQ(ierr);
535           } else {
536             ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,nc,work+cum,&S_FE);CHKERRQ(ierr);
537             cum  = nc*subset_size;
538             ierr = MatGetSubMatrix(S_all,is_E,is_F,MAT_REUSE_MATRIX,&S_EF);CHKERRQ(ierr);
539           }
540           ierr = ISDestroy(&is_F);CHKERRQ(ierr);
541           if (chop) {
542             ierr = MatChop(S_FF,eps);CHKERRQ(ierr);
543             ierr = MatConvert(S_FF,MATAIJ,MAT_REUSE_MATRIX,&S_FF);CHKERRQ(ierr);
544           }
545           ierr = MatCreateSchurComplement(S_FF,S_FF,S_FE,S_EF,sub_schurs->S_Ej[i],&Stilda_impl);CHKERRQ(ierr);
546           ierr = MatDestroy(&S_FF);CHKERRQ(ierr);
547           ierr = MatDestroy(&S_FE);CHKERRQ(ierr);
548           ierr = MatDestroy(&S_EF);CHKERRQ(ierr);
549           ierr = MatSchurComplementGetKSP(Stilda_impl,&ksp);CHKERRQ(ierr);
550           ierr = KSPSetType(ksp,KSPPREONLY);CHKERRQ(ierr);
551           ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr);
552           if (sub_schurs->is_hermitian && sub_schurs->is_posdef) {
553             ierr = PCSetType(pc,PCCHOLESKY);CHKERRQ(ierr);
554           } else {
555             ierr = PCSetType(pc,PCLU);CHKERRQ(ierr);
556           }
557           if (!chop) {
558             ierr = PCFactorSetUseInPlace(pc,PETSC_TRUE);CHKERRQ(ierr);
559           } else {
560             ierr = PCFactorSetMatSolverPackage(pc,MATSOLVERMUMPS);CHKERRQ(ierr);
561           }
562           ierr = KSPSetUp(ksp);CHKERRQ(ierr);
563           ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work+cum,&Stilda);CHKERRQ(ierr);
564           ierr = PCBDDCComputeExplicitSchur(Stilda_impl,sub_schurs->is_hermitian,MAT_REUSE_MATRIX,&Stilda);CHKERRQ(ierr);
565           ierr = MatDestroy(&Stilda_impl);CHKERRQ(ierr);
566         }
567 /*
568         PetscViewerPushFormat(PETSC_VIEWER_STDOUT_SELF,PETSC_VIEWER_ASCII_MATLAB);
569         ierr = MatView(Stilda,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
570         PetscViewerPopFormat(PETSC_VIEWER_STDOUT_SELF);
571 */
572         ierr = MatDenseGetArray(Stilda,&vals);CHKERRQ(ierr);
573         if (deluxe) { /* when using deluxe scaling, we need (S_1^-1+S_2^-1)^-1 */
574           PetscScalar *vals2;
575 
576           ierr = MatDenseGetArray(sub_schurs->S_Ej[i],&vals2);CHKERRQ(ierr);
577           /* need to be optimized (cholesky) */
578           ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr);
579           ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
580           if (invert_Stildas) { /* Stildas can be singular */
581             PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,vals,&B_N,pivots,&B_ierr));
582             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
583             PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,vals,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
584             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
585           }
586           PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,vals2,&B_N,pivots,&B_ierr));
587           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
588           PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,vals2,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
589           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
590           ierr = PetscFPTrapPop();CHKERRQ(ierr);
591           ierr = MatSetValues(S_Ej_inv_all,subset_size,dummy_idx,subset_size,dummy_idx,vals2,INSERT_VALUES);CHKERRQ(ierr);
592           ierr = MatDenseRestoreArray(sub_schurs->S_Ej[i],&vals2);CHKERRQ(ierr);
593         }
594         ierr = MatSetValues(S_Ej_tilda_all,subset_size,dummy_idx,subset_size,dummy_idx,vals,INSERT_VALUES);CHKERRQ(ierr);
595         ierr = MatDenseRestoreArray(Stilda,&vals);CHKERRQ(ierr);
596         ierr = MatDestroy(&Stilda);CHKERRQ(ierr);
597       }
598       ierr = ISDestroy(&is_E);CHKERRQ(ierr);
599       local_size += subset_size;
600     }
601     ierr = PetscFree2(dummy_idx,work);CHKERRQ(ierr);
602   }
603   ierr = PetscFree(nnz);CHKERRQ(ierr);
604   ierr = MatDestroy(&S_all);CHKERRQ(ierr);
605   ierr = MatAssemblyBegin(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
606   ierr = MatAssemblyEnd(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
607   if (S_Ej_tilda_all) {
608     ierr = MatAssemblyBegin(S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
609     ierr = MatAssemblyEnd(S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
610   }
611   if (S_Ej_inv_all) {
612     ierr = MatAssemblyBegin(S_Ej_inv_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
613     ierr = MatAssemblyEnd(S_Ej_inv_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
614   }
615 
616   /* Global matrix of all assembled Schur on subsets */
617   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);
618   ierr = ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)sub_schurs->l2gmap),1,local_size,all_local_idx_G,PETSC_COPY_VALUES,&l2gmap_subsets);CHKERRQ(ierr);
619   ierr = MatCreateIS(PetscObjectComm((PetscObject)sub_schurs->l2gmap),1,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size,l2gmap_subsets,&work_mat);CHKERRQ(ierr);
620   ierr = ISLocalToGlobalMappingDestroy(&l2gmap_subsets);CHKERRQ(ierr);
621   ierr = MatISSetLocalMat(work_mat,sub_schurs->S_Ej_all);CHKERRQ(ierr);
622   ierr = MatCreate(PetscObjectComm((PetscObject)work_mat),&global_schur_subsets);CHKERRQ(ierr);
623   ierr = MatSetSizes(global_schur_subsets,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size);CHKERRQ(ierr);
624   ierr = MatSetType(global_schur_subsets,MATMPIAIJ);CHKERRQ(ierr);
625   ierr = MatISSetMPIXAIJPreallocation_Private(work_mat,global_schur_subsets,PETSC_TRUE);CHKERRQ(ierr);
626   ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
627   /* Get local part of (\sum_j S_Ej) */
628   ierr = PetscFree(all_local_idx_N);CHKERRQ(ierr);
629   ierr = ISCreateGeneral(PetscObjectComm((PetscObject)sub_schurs->l2gmap),local_size,all_local_idx_G,PETSC_OWN_POINTER,&temp_is);CHKERRQ(ierr);
630   if (sub_schurs->sum_S_Ej_all) {
631     ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,MAT_REUSE_MATRIX,&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
632   } else {
633     ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,MAT_INITIAL_MATRIX,&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
634   }
635 
636   if (S_Ej_tilda_all) {
637     ierr = MatISSetLocalMat(work_mat,S_Ej_tilda_all);CHKERRQ(ierr);
638     ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
639     if (sub_schurs->sum_S_Ej_tilda_all) {
640       ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,MAT_REUSE_MATRIX,&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
641     } else {
642       ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,MAT_INITIAL_MATRIX,&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
643     }
644   }
645   if (S_Ej_inv_all) {
646     PetscInt    cum;
647     PetscScalar *array,*array2;
648     ierr = MatISSetLocalMat(work_mat,S_Ej_inv_all);CHKERRQ(ierr);
649     ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
650     if (sub_schurs->sum_S_Ej_inv_all) {
651       ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,MAT_REUSE_MATRIX,&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
652     } else {
653       ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,MAT_INITIAL_MATRIX,&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
654     }
655     /* invert blocks of sum_S_Ej_inv_all */
656     ierr = MatSeqAIJGetArray(sub_schurs->sum_S_Ej_inv_all,&array);CHKERRQ(ierr);
657     ierr = MatSeqAIJGetArray(sub_schurs->sum_S_Ej_tilda_all,&array2);CHKERRQ(ierr);
658     cum = 0;
659     for (i=0;i<sub_schurs->n_subs;i++) {
660       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
661       if (PetscBTLookup(sub_schurs->computed_Stilda_subs,i)) {
662         /* need to be optimized (cholesky) */
663         ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr);
664         ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
665         PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,array+cum,&B_N,pivots,&B_ierr));
666         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
667         PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,array+cum,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
668         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
669         if (invert_Stildas) { /* Stildas can be singular */
670           PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,array2+cum,&B_N,pivots,&B_ierr));
671           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
672           PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,array2+cum,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
673           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
674         }
675         ierr = PetscFPTrapPop();CHKERRQ(ierr);
676       }
677       cum += subset_size*subset_size;
678     }
679     ierr = MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_inv_all,&array);CHKERRQ(ierr);
680     ierr = MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_tilda_all,&array2);CHKERRQ(ierr);
681   }
682 
683   ierr = PetscFree2(Bwork,pivots);CHKERRQ(ierr);
684   ierr = MatDestroy(&global_schur_subsets);CHKERRQ(ierr);
685   ierr = MatDestroy(&S_Ej_tilda_all);CHKERRQ(ierr);
686   ierr = MatDestroy(&S_Ej_inv_all);CHKERRQ(ierr);
687   ierr = MatDestroy(&work_mat);CHKERRQ(ierr);
688   ierr = ISDestroy(&temp_is);CHKERRQ(ierr);
689   PetscFunctionReturn(0);
690 }
691 
692 #undef __FUNCT__
693 #define __FUNCT__ "PCBDDCSubSchursInit"
694 PetscErrorCode PCBDDCSubSchursInit(PCBDDCSubSchurs sub_schurs, Mat A, Mat S, IS is_I, IS is_B, PCBDDCGraph graph, ISLocalToGlobalMapping BtoNmap, PetscInt seqthreshold)
695 {
696   IS                  *faces,*edges,*all_cc,vertices;
697   PetscInt            *index_sequential,*index_parallel;
698   PetscInt            *auxlocal_sequential,*auxlocal_parallel;
699   PetscInt            *auxglobal_sequential,*auxglobal_parallel;
700   PetscInt            *auxmapping;
701   PetscInt            i,max_subset_size;
702   PetscInt            n_sequential_problems,n_local_sequential_problems,n_parallel_problems,n_local_parallel_problems;
703   PetscInt            n_faces,n_edges,n_all_cc;
704   PetscBool           is_sorted;
705   PetscErrorCode  ierr;
706 
707   PetscFunctionBegin;
708   ierr = ISSorted(is_I,&is_sorted);CHKERRQ(ierr);
709   if (!is_sorted) {
710     SETERRQ(PetscObjectComm((PetscObject)is_I),PETSC_ERR_PLIB,"IS for I dofs should be shorted");
711   }
712   ierr = ISSorted(is_B,&is_sorted);CHKERRQ(ierr);
713   if (!is_sorted) {
714     SETERRQ(PetscObjectComm((PetscObject)is_B),PETSC_ERR_PLIB,"IS for B dofs should be shorted");
715   }
716 
717   /* reset any previous data */
718   ierr = PCBDDCSubSchursReset(sub_schurs);CHKERRQ(ierr);
719 
720   /* get index sets for faces and edges */
721   ierr = PCBDDCGraphGetCandidatesIS(graph,&n_faces,&faces,&n_edges,&edges,&vertices);CHKERRQ(ierr);
722   n_all_cc = n_faces+n_edges;
723   ierr = PetscBTCreate(n_all_cc,&sub_schurs->is_edge);CHKERRQ(ierr);
724   ierr = PetscBTCreate(n_all_cc,&sub_schurs->computed_Stilda_subs);CHKERRQ(ierr);
725   ierr = PetscMalloc1(n_all_cc,&all_cc);CHKERRQ(ierr);
726   for (i=0;i<n_faces;i++) {
727     all_cc[i] = faces[i];
728   }
729   for (i=0;i<n_edges;i++) {
730     all_cc[n_faces+i] = edges[i];
731     ierr = PetscBTSet(sub_schurs->is_edge,n_faces+i);CHKERRQ(ierr);
732   }
733   ierr = PetscFree(faces);CHKERRQ(ierr);
734   ierr = PetscFree(edges);CHKERRQ(ierr);
735 
736   /* map interface's subsets */
737   max_subset_size = 0;
738   for (i=0;i<n_all_cc;i++) {
739     PetscInt subset_size;
740     ierr = ISGetLocalSize(all_cc[i],&subset_size);CHKERRQ(ierr);
741     max_subset_size = PetscMax(max_subset_size,subset_size);
742   }
743   ierr = PetscMalloc1(max_subset_size,&auxmapping);CHKERRQ(ierr);
744   ierr = PetscMalloc2(graph->ncc,&auxlocal_sequential,
745                       graph->ncc,&auxlocal_parallel);CHKERRQ(ierr);
746   ierr = PetscMalloc2(graph->ncc,&index_sequential,
747                       graph->ncc,&index_parallel);CHKERRQ(ierr);
748 
749   /* if threshold is negative uses all sequential problems (possibly using MUMPS) */
750   sub_schurs->use_mumps = PETSC_FALSE;
751   if (seqthreshold < 0) {
752     seqthreshold = max_subset_size;
753 #if defined(PETSC_HAVE_MUMPS)
754     sub_schurs->use_mumps = !!A;
755 #endif
756   }
757 
758 
759   /* determine which problem has to be solved in parallel or sequentially */
760   n_local_sequential_problems = 0;
761   n_local_parallel_problems = 0;
762   for (i=0;i<n_all_cc;i++) {
763     PetscInt       subset_size,j,min_loc = 0;
764     const PetscInt *idxs;
765 
766     ierr = ISGetLocalSize(all_cc[i],&subset_size);CHKERRQ(ierr);
767     ierr = ISGetIndices(all_cc[i],&idxs);CHKERRQ(ierr);
768     ierr = ISLocalToGlobalMappingApply(graph->l2gmap,subset_size,idxs,auxmapping);CHKERRQ(ierr);
769     for (j=1;j<subset_size;j++) {
770       if (auxmapping[j]<auxmapping[min_loc]) {
771         min_loc = j;
772       }
773     }
774     if (subset_size > seqthreshold) {
775       index_parallel[n_local_parallel_problems] = i;
776       auxlocal_parallel[n_local_parallel_problems] = idxs[min_loc];
777       n_local_parallel_problems++;
778     } else {
779       index_sequential[n_local_sequential_problems] = i;
780       auxlocal_sequential[n_local_sequential_problems] = idxs[min_loc];
781       n_local_sequential_problems++;
782     }
783     ierr = ISRestoreIndices(all_cc[i],&idxs);CHKERRQ(ierr);
784   }
785 
786   /* Number parallel problems */
787   auxglobal_parallel = 0;
788   ierr = PCBDDCSubsetNumbering(PetscObjectComm((PetscObject)graph->l2gmap),graph->l2gmap,n_local_parallel_problems,auxlocal_parallel,PETSC_NULL,&n_parallel_problems,&auxglobal_parallel);CHKERRQ(ierr);
789 
790   /* Number sequential problems */
791   auxglobal_sequential = 0;
792   ierr = PCBDDCSubsetNumbering(PetscObjectComm((PetscObject)graph->l2gmap),graph->l2gmap,n_local_sequential_problems,auxlocal_sequential,PETSC_NULL,&n_sequential_problems,&auxglobal_sequential);CHKERRQ(ierr);
793 
794   /* update info in sub_schurs */
795   if (A) {
796     PetscBool isseqaij;
797 
798     ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
799     if (isseqaij) {
800       ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
801       sub_schurs->A = A;
802     } else { /* SeqBAIJ matrices does not support symmetry checking, SeqSBAIJ does not support MatPermute */
803       ierr = MatConvert(A,MATSEQAIJ,MAT_INITIAL_MATRIX,&sub_schurs->A);CHKERRQ(ierr);
804     }
805   }
806   ierr = PetscObjectReference((PetscObject)S);CHKERRQ(ierr);
807   sub_schurs->S = S;
808   ierr = PetscObjectReference((PetscObject)is_I);CHKERRQ(ierr);
809   sub_schurs->is_I = is_I;
810   ierr = PetscObjectReference((PetscObject)is_B);CHKERRQ(ierr);
811   sub_schurs->is_B = is_B;
812   ierr = PetscObjectReference((PetscObject)graph->l2gmap);CHKERRQ(ierr);
813   sub_schurs->l2gmap = graph->l2gmap;
814   ierr = PetscObjectReference((PetscObject)BtoNmap);CHKERRQ(ierr);
815   sub_schurs->BtoNmap = BtoNmap;
816   sub_schurs->n_subs_seq = n_local_sequential_problems;
817   sub_schurs->n_subs_par = n_local_parallel_problems;
818   sub_schurs->n_subs_seq_g = n_sequential_problems;
819   sub_schurs->n_subs_par_g = n_parallel_problems;
820   sub_schurs->n_subs = sub_schurs->n_subs_seq + sub_schurs->n_subs_par;
821   sub_schurs->is_subs = all_cc;
822   if (!sub_schurs->use_mumps) { /* for adaptive selection */
823     for (i=0;i<sub_schurs->n_subs;i++) {
824       ierr = ISSort(sub_schurs->is_subs[i]);CHKERRQ(ierr);
825     }
826   }
827   sub_schurs->is_Ej_com = vertices;
828   sub_schurs->index_sequential = index_sequential;
829   sub_schurs->index_parallel = index_parallel;
830   sub_schurs->auxglobal_sequential = auxglobal_sequential;
831   sub_schurs->auxglobal_parallel = auxglobal_parallel;
832 
833 
834   /* allocate space for schur complements */
835   ierr = PetscCalloc1(sub_schurs->n_subs,&sub_schurs->S_Ej);CHKERRQ(ierr);
836   sub_schurs->S_Ej_all = NULL;
837   sub_schurs->sum_S_Ej_all = NULL;
838   sub_schurs->sum_S_Ej_inv_all = NULL;
839   sub_schurs->sum_S_Ej_tilda_all = NULL;
840   sub_schurs->is_Ej_all = NULL;
841 
842   /* free workspace */
843   ierr = PetscFree(auxmapping);CHKERRQ(ierr);
844   ierr = PetscFree2(auxlocal_sequential,auxlocal_parallel);CHKERRQ(ierr);
845 
846   PetscFunctionReturn(0);
847 }
848 
849 #undef __FUNCT__
850 #define __FUNCT__ "PCBDDCSubSchursCreate"
851 PetscErrorCode PCBDDCSubSchursCreate(PCBDDCSubSchurs *sub_schurs)
852 {
853   PCBDDCSubSchurs schurs_ctx;
854   PetscErrorCode  ierr;
855 
856   PetscFunctionBegin;
857   ierr = PetscNew(&schurs_ctx);CHKERRQ(ierr);
858   schurs_ctx->n_subs = 0;
859   *sub_schurs = schurs_ctx;
860   PetscFunctionReturn(0);
861 }
862 
863 #undef __FUNCT__
864 #define __FUNCT__ "PCBDDCSubSchursDestroy"
865 PetscErrorCode PCBDDCSubSchursDestroy(PCBDDCSubSchurs *sub_schurs)
866 {
867   PetscErrorCode ierr;
868 
869   PetscFunctionBegin;
870   ierr = PCBDDCSubSchursReset(*sub_schurs);CHKERRQ(ierr);
871   ierr = PetscFree(*sub_schurs);CHKERRQ(ierr);
872   PetscFunctionReturn(0);
873 }
874 
875 #undef __FUNCT__
876 #define __FUNCT__ "PCBDDCSubSchursReset"
877 PetscErrorCode PCBDDCSubSchursReset(PCBDDCSubSchurs sub_schurs)
878 {
879   PetscInt       i;
880   PetscErrorCode ierr;
881 
882   PetscFunctionBegin;
883   ierr = MatDestroy(&sub_schurs->A);CHKERRQ(ierr);
884   ierr = MatDestroy(&sub_schurs->S);CHKERRQ(ierr);
885   ierr = ISDestroy(&sub_schurs->is_I);CHKERRQ(ierr);
886   ierr = ISDestroy(&sub_schurs->is_B);CHKERRQ(ierr);
887   ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->l2gmap);CHKERRQ(ierr);
888   ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->BtoNmap);CHKERRQ(ierr);
889   ierr = MatDestroy(&sub_schurs->S_Ej_all);CHKERRQ(ierr);
890   ierr = MatDestroy(&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
891   ierr = MatDestroy(&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
892   ierr = MatDestroy(&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
893   ierr = ISDestroy(&sub_schurs->is_Ej_all);CHKERRQ(ierr);
894   ierr = ISDestroy(&sub_schurs->is_Ej_com);CHKERRQ(ierr);
895   ierr = PetscBTDestroy(&sub_schurs->is_edge);CHKERRQ(ierr);
896   ierr = PetscBTDestroy(&sub_schurs->computed_Stilda_subs);CHKERRQ(ierr);
897   for (i=0;i<sub_schurs->n_subs;i++) {
898     ierr = ISDestroy(&sub_schurs->is_subs[i]);CHKERRQ(ierr);
899     ierr = MatDestroy(&sub_schurs->S_Ej[i]);CHKERRQ(ierr);
900   }
901   ierr = ISDestroy(&sub_schurs->is_I_layer);CHKERRQ(ierr);
902   if (sub_schurs->n_subs) {
903     ierr = PetscFree(sub_schurs->is_subs);CHKERRQ(ierr);
904     ierr = PetscFree(sub_schurs->S_Ej);CHKERRQ(ierr);
905   }
906   ierr = PetscFree2(sub_schurs->index_sequential,sub_schurs->index_parallel);CHKERRQ(ierr);
907   ierr = PetscFree(sub_schurs->auxglobal_sequential);CHKERRQ(ierr);
908   ierr = PetscFree(sub_schurs->auxglobal_parallel);CHKERRQ(ierr);
909   sub_schurs->n_subs = 0;
910   PetscFunctionReturn(0);
911 }
912 
913 #undef __FUNCT__
914 #define __FUNCT__ "PCBDDCAdjGetNextLayer_Private"
915 PETSC_STATIC_INLINE PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt* queue_tip,PetscInt n_prev,PetscBT touched,PetscInt* xadj,PetscInt* adjncy,PetscInt* n_added)
916 {
917   PetscInt       i,j,n;
918   PetscErrorCode ierr;
919 
920   PetscFunctionBegin;
921   n = 0;
922   for (i=-n_prev;i<0;i++) {
923     PetscInt start_dof = queue_tip[i];
924     for (j=xadj[start_dof];j<xadj[start_dof+1];j++) {
925       PetscInt dof = adjncy[j];
926       if (!PetscBTLookup(touched,dof)) {
927         ierr = PetscBTSet(touched,dof);CHKERRQ(ierr);
928         queue_tip[n] = dof;
929         n++;
930       }
931     }
932   }
933   *n_added = n;
934   PetscFunctionReturn(0);
935 }
936