xref: /petsc/src/ksp/pc/impls/bddc/bddcschurs.c (revision 2b9730978a4291cabe2a55e15026fd06766dc70e)
1 #include <../src/ksp/pc/impls/bddc/bddc.h>
2 #include <../src/ksp/pc/impls/bddc/bddcprivate.h>
3 
4 PETSC_STATIC_INLINE PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt*,PetscInt,PetscBT,PetscInt*,PetscInt*,PetscInt*);
5 static PetscErrorCode PCBDDCComputeExplicitSchur(Mat M, Mat *S);
6 
7 #undef __FUNCT__
8 #define __FUNCT__ "PCBDDCComputeExplicitSchur"
9 static PetscErrorCode PCBDDCComputeExplicitSchur(Mat M, Mat *S)
10 {
11   Mat            B, C, D, Bd, Cd, AinvBd;
12   KSP            ksp;
13   PC             pc;
14   PetscBool      isLU, isILU, isCHOL, Bdense, Cdense;
15   PetscReal      fill = 2.0;
16   PetscMPIInt    size;
17   PetscErrorCode ierr;
18 
19   PetscFunctionBegin;
20   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)M),&size);CHKERRQ(ierr);
21   if (size != 1) {
22     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not for parallel matrices");
23   }
24   ierr = MatSchurComplementGetSubMatrices(M, NULL, NULL, &B, &C, &D);CHKERRQ(ierr);
25   ierr = MatSchurComplementGetKSP(M, &ksp);CHKERRQ(ierr);
26   ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr);
27   ierr = PetscObjectTypeCompare((PetscObject) pc, PCLU, &isLU);CHKERRQ(ierr);
28   ierr = PetscObjectTypeCompare((PetscObject) pc, PCILU, &isILU);CHKERRQ(ierr);
29   ierr = PetscObjectTypeCompare((PetscObject) pc, PCCHOLESKY, &isCHOL);CHKERRQ(ierr);
30   ierr = PetscObjectTypeCompare((PetscObject) B, MATSEQDENSE, &Bdense);CHKERRQ(ierr);
31   ierr = PetscObjectTypeCompare((PetscObject) C, MATSEQDENSE, &Cdense);CHKERRQ(ierr);
32   if (!Bdense) {
33     ierr = MatConvert(B, MATSEQDENSE, MAT_INITIAL_MATRIX, &Bd);CHKERRQ(ierr);
34   } else {
35     Bd = B;
36   }
37 
38   if (isLU || isILU || isCHOL) {
39     Mat fact;
40 
41     ierr = KSPSetUp(ksp);CHKERRQ(ierr);
42     ierr = PCFactorGetMatrix(pc, &fact);CHKERRQ(ierr);
43     ierr = MatDuplicate(Bd, MAT_DO_NOT_COPY_VALUES, &AinvBd);CHKERRQ(ierr);
44     ierr = MatMatSolve(fact, Bd, AinvBd);CHKERRQ(ierr);
45   } else {
46     Mat Ainvd;
47 
48     ierr = PCComputeExplicitOperator(pc, &Ainvd);CHKERRQ(ierr);
49     ierr = MatMatMult(Ainvd, Bd, MAT_INITIAL_MATRIX, fill, &AinvBd);CHKERRQ(ierr);
50     ierr = MatDestroy(&Ainvd);CHKERRQ(ierr);
51   }
52   if (!Bdense) {
53     ierr = MatDestroy(&Bd);CHKERRQ(ierr);
54   }
55   if (!Cdense) {
56     ierr = MatConvert(C, MATSEQDENSE, MAT_INITIAL_MATRIX, &Cd);CHKERRQ(ierr);
57   } else {
58     Cd = C;
59   }
60 
61   ierr = MatMatMult(Cd, AinvBd, MAT_INITIAL_MATRIX, fill, S);CHKERRQ(ierr);
62   ierr = MatDestroy(&AinvBd);CHKERRQ(ierr);
63   if (!Cdense) {
64     ierr = MatDestroy(&Cd);CHKERRQ(ierr);
65   }
66 
67   if (D) {
68     Mat       Dd;
69     PetscBool Ddense;
70 
71     ierr = PetscObjectTypeCompare((PetscObject)D,MATSEQDENSE,&Ddense);CHKERRQ(ierr);
72     if (!Ddense) {
73       ierr = MatConvert(D, MATSEQDENSE, MAT_INITIAL_MATRIX, &Dd);CHKERRQ(ierr);
74     } else {
75       Dd = D;
76     }
77     ierr = MatAYPX(*S,-1.0,Dd,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
78     if (!Ddense) {
79       ierr = MatDestroy(&Dd);CHKERRQ(ierr);
80     }
81   } else {
82     ierr = MatScale(*S,-1.0);CHKERRQ(ierr);
83   }
84   PetscFunctionReturn(0);
85 }
86 
87 #undef __FUNCT__
88 #define __FUNCT__ "PCBDDCSubSchursSetUp"
89 PetscErrorCode PCBDDCSubSchursSetUp(PCBDDCSubSchurs sub_schurs, PetscInt xadj[], PetscInt adjncy[], PetscInt nlayers)
90 {
91   Mat                    A_II,A_IB,A_BI,A_BB;
92   Mat                    AE_II,*AE_IE,*AE_EI,*AE_EE;
93   Mat                    S_all,global_schur_subsets,work_mat;
94   ISLocalToGlobalMapping l2gmap_subsets;
95   IS                     is_I,*is_subset_B,temp_is;
96   PetscInt               *nnz,*all_local_idx_G,*all_local_idx_B,*all_local_idx_N;
97   PetscInt               i,subset_size,max_subset_size;
98   PetscInt               extra,local_size,global_size;
99   PetscErrorCode         ierr;
100 
101   PetscFunctionBegin;
102 
103   /* allocate space for schur complements */
104   ierr = PetscMalloc5(sub_schurs->n_subs,&sub_schurs->is_AEj_B,
105                       sub_schurs->n_subs,&sub_schurs->S_Ej,
106                       sub_schurs->n_subs,&sub_schurs->S_Ej_tilda,
107                       sub_schurs->n_subs,&sub_schurs->work1,
108                       sub_schurs->n_subs,&sub_schurs->work2);CHKERRQ(ierr);
109 
110   /* get Schur complement matrices */
111   if (!sub_schurs->use_mumps) {
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   }
118 
119   /* determine interior problems */
120   if (nlayers >= 0 && xadj != NULL && adjncy != NULL) { /* Interior problems can be different from the original one */
121     PetscBT                touched;
122     const PetscInt*        idx_B;
123     PetscInt               n_I,n_B,n_local_dofs,n_prev_added,j,layer,*local_numbering;
124 
125     /* get sizes */
126     ierr = ISGetLocalSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr);
127     ierr = ISGetLocalSize(sub_schurs->is_B,&n_B);CHKERRQ(ierr);
128 
129     ierr = PetscMalloc1(n_I+n_B,&local_numbering);CHKERRQ(ierr);
130     ierr = PetscBTCreate(n_I+n_B,&touched);CHKERRQ(ierr);
131     ierr = PetscBTMemzero(n_I+n_B,touched);CHKERRQ(ierr);
132 
133     /* all boundary dofs must be skipped when adding layers */
134     ierr = ISGetIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr);
135     for (j=0;j<n_B;j++) {
136       ierr = PetscBTSet(touched,idx_B[j]);CHKERRQ(ierr);
137     }
138     ierr = PetscMemcpy(local_numbering,idx_B,n_B*sizeof(PetscInt));CHKERRQ(ierr);
139     ierr = ISRestoreIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr);
140 
141     /* add prescribed number of layers of dofs */
142     n_local_dofs = n_B;
143     n_prev_added = n_B;
144     for (layer=0;layer<nlayers;layer++) {
145       PetscInt n_added;
146       if (n_local_dofs == n_I+n_B) break;
147       if (n_local_dofs > n_I+n_B) {
148         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);
149       }
150       ierr = PCBDDCAdjGetNextLayer_Private(local_numbering+n_local_dofs,n_prev_added,touched,xadj,adjncy,&n_added);CHKERRQ(ierr);
151       n_prev_added = n_added;
152       n_local_dofs += n_added;
153       if (!n_added) break;
154     }
155     ierr = PetscBTDestroy(&touched);CHKERRQ(ierr);
156 
157     /* IS for I layer dofs in original numbering */
158     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);
159     ierr = PetscFree(local_numbering);CHKERRQ(ierr);
160     ierr = ISSort(sub_schurs->is_I_layer);CHKERRQ(ierr);
161     /* IS for I layer dofs in I numbering */
162     if (!sub_schurs->use_mumps) {
163       ISLocalToGlobalMapping ItoNmap;
164       ierr = ISLocalToGlobalMappingCreateIS(sub_schurs->is_I,&ItoNmap);CHKERRQ(ierr);
165       ierr = ISGlobalToLocalMappingApplyIS(ItoNmap,IS_GTOLM_DROP,sub_schurs->is_I_layer,&is_I);CHKERRQ(ierr);
166       ierr = ISLocalToGlobalMappingDestroy(&ItoNmap);CHKERRQ(ierr);
167 
168       /* II block */
169       ierr = MatGetSubMatrix(A_II,is_I,is_I,MAT_INITIAL_MATRIX,&AE_II);CHKERRQ(ierr);
170     }
171   } else {
172     PetscInt n_I;
173 
174     /* IS for I dofs in original numbering */
175     ierr = PetscObjectReference((PetscObject)sub_schurs->is_I);CHKERRQ(ierr);
176     sub_schurs->is_I_layer = sub_schurs->is_I;
177 
178     /* IS for I dofs in I numbering (strided 1) */
179     if (!sub_schurs->use_mumps) {
180       ierr = ISGetSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr);
181       ierr = ISCreateStride(PetscObjectComm((PetscObject)sub_schurs->is_I),n_I,0,1,&is_I);CHKERRQ(ierr);
182 
183       /* II block is the same */
184       ierr = PetscObjectReference((PetscObject)A_II);CHKERRQ(ierr);
185       AE_II = A_II;
186     }
187   }
188 
189   for (i=0;i<sub_schurs->n_subs;i++) {
190     ierr = ISDuplicate(sub_schurs->is_subs[i],&sub_schurs->is_AEj_B[i]);CHKERRQ(ierr);
191     ierr = ISSort(sub_schurs->is_AEj_B[i]);CHKERRQ(ierr);
192   }
193 
194   /* Get info on subset sizes and sum of all subsets sizes */
195   max_subset_size = 0;
196   local_size = 0;
197   for (i=0;i<sub_schurs->n_subs_seq;i++) {
198     PetscInt j = sub_schurs->index_sequential[i];
199     ierr = ISGetLocalSize(sub_schurs->is_AEj_B[j],&subset_size);CHKERRQ(ierr);
200     max_subset_size = PetscMax(subset_size,max_subset_size);
201     local_size += subset_size;
202   }
203 
204   /* Work arrays for local indices */
205   ierr = PetscMalloc1(local_size,&all_local_idx_B);CHKERRQ(ierr);
206   extra = 0;
207   if (sub_schurs->use_mumps) {
208     ierr = ISGetLocalSize(sub_schurs->is_I_layer,&extra);CHKERRQ(ierr);
209   }
210   ierr = PetscMalloc1(local_size+extra,&all_local_idx_N);CHKERRQ(ierr);
211   if (extra) {
212     const PetscInt *idxs;
213     ierr = ISGetIndices(sub_schurs->is_I_layer,&idxs);CHKERRQ(ierr);
214     ierr = PetscMemcpy(all_local_idx_N,idxs,extra*sizeof(PetscInt));CHKERRQ(ierr);
215     ierr = ISRestoreIndices(sub_schurs->is_I_layer,&idxs);CHKERRQ(ierr);
216   }
217   ierr = PetscMalloc1(local_size,&nnz);CHKERRQ(ierr);
218 
219   /* Get local indices in local numbering */
220   local_size = 0;
221   for (i=0;i<sub_schurs->n_subs_seq;i++) {
222     PetscInt j;
223     const    PetscInt *idxs;
224 
225     PetscInt local_problem_index = sub_schurs->index_sequential[i];
226     ierr = ISGetLocalSize(sub_schurs->is_AEj_B[local_problem_index],&subset_size);CHKERRQ(ierr);
227     ierr = ISGetIndices(sub_schurs->is_AEj_B[local_problem_index],&idxs);CHKERRQ(ierr);
228     /* subset indices in local numbering */
229     ierr = PetscMemcpy(all_local_idx_N+local_size+extra,idxs,subset_size*sizeof(PetscInt));CHKERRQ(ierr);
230     ierr = ISRestoreIndices(sub_schurs->is_AEj_B[local_problem_index],&idxs);CHKERRQ(ierr);
231     for (j=0;j<subset_size;j++) nnz[local_size+j] = subset_size;
232     local_size += subset_size;
233   }
234 
235   S_all = NULL;
236   if (!sub_schurs->use_mumps) {
237     /* subsets in original and boundary numbering */
238     for (i=0;i<sub_schurs->n_subs;i++) {
239       ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,sub_schurs->is_AEj_B[i],&is_subset_B[i]);CHKERRQ(ierr);
240     }
241 
242     /* EE block */
243     for (i=0;i<sub_schurs->n_subs;i++) {
244       ierr = MatGetSubMatrix(A_BB,is_subset_B[i],is_subset_B[i],MAT_INITIAL_MATRIX,&AE_EE[i]);CHKERRQ(ierr);
245     }
246     /* IE block */
247     for (i=0;i<sub_schurs->n_subs;i++) {
248       ierr = MatGetSubMatrix(A_IB,is_I,is_subset_B[i],MAT_INITIAL_MATRIX,&AE_IE[i]);CHKERRQ(ierr);
249     }
250     /* EI block */
251     for (i=0;i<sub_schurs->n_subs;i++) {
252       ierr = MatGetSubMatrix(A_BI,is_subset_B[i],is_I,MAT_INITIAL_MATRIX,&AE_EI[i]);CHKERRQ(ierr);
253     }
254 
255     /* setup Schur complements on subset */
256     for (i=0;i<sub_schurs->n_subs;i++) {
257       ierr = MatCreateSchurComplement(AE_II,AE_II,AE_IE[i],AE_EI[i],AE_EE[i],&sub_schurs->S_Ej[i]);CHKERRQ(ierr);
258       if (AE_II == A_II) { /* we can reuse the same ksp */
259         KSP ksp;
260         ierr = MatSchurComplementGetKSP(sub_schurs->S,&ksp);CHKERRQ(ierr);
261         ierr = MatSchurComplementSetKSP(sub_schurs->S_Ej[i],ksp);CHKERRQ(ierr);
262       } else { /* build new ksp object which inherits ksp and pc types from the original one */
263         KSP      origksp,schurksp;
264         PC       origpc,schurpc;
265         KSPType  ksp_type;
266         PCType   pc_type;
267         PetscInt n_internal;
268 
269         ierr = MatSchurComplementGetKSP(sub_schurs->S,&origksp);CHKERRQ(ierr);
270         ierr = MatSchurComplementGetKSP(sub_schurs->S_Ej[i],&schurksp);CHKERRQ(ierr);
271         ierr = KSPGetType(origksp,&ksp_type);CHKERRQ(ierr);
272         ierr = KSPSetType(schurksp,ksp_type);CHKERRQ(ierr);
273         ierr = KSPGetPC(schurksp,&schurpc);CHKERRQ(ierr);
274         ierr = KSPGetPC(origksp,&origpc);CHKERRQ(ierr);
275         ierr = PCGetType(origpc,&pc_type);CHKERRQ(ierr);
276         ierr = PCSetType(schurpc,pc_type);CHKERRQ(ierr);
277         ierr = ISGetSize(is_I,&n_internal);CHKERRQ(ierr);
278         if (n_internal) { /* UMFPACK gives error with 0 sized problems */
279           MatSolverPackage solver=NULL;
280           ierr = PCFactorGetMatSolverPackage(origpc,(const MatSolverPackage*)&solver);CHKERRQ(ierr);
281           if (solver) {
282             ierr = PCFactorSetMatSolverPackage(schurpc,solver);CHKERRQ(ierr);
283           }
284         }
285         ierr = KSPSetUp(schurksp);CHKERRQ(ierr);
286       }
287     }
288     /* free */
289     ierr = ISDestroy(&is_I);CHKERRQ(ierr);
290     ierr = MatDestroy(&AE_II);CHKERRQ(ierr);
291     for (i=0;i<sub_schurs->n_subs;i++) {
292       ierr = MatDestroy(&AE_EE[i]);CHKERRQ(ierr);
293       ierr = MatDestroy(&AE_IE[i]);CHKERRQ(ierr);
294       ierr = MatDestroy(&AE_EI[i]);CHKERRQ(ierr);
295       ierr = ISDestroy(&is_subset_B[i]);CHKERRQ(ierr);
296     }
297     ierr = PetscFree4(is_subset_B,AE_IE,AE_EI,AE_EE);CHKERRQ(ierr);
298 #if defined(PETSC_HAVE_MUMPS)
299   } else {
300     Mat           A,F;
301     IS            is_A_all;
302     PetscBool     is_symmetric;
303     PetscInt      *idxs_schur,n_I;
304 
305     /* get working mat */
306     ierr = ISGetLocalSize(sub_schurs->is_I_layer,&n_I);CHKERRQ(ierr);
307     ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size+n_I,all_local_idx_N,PETSC_COPY_VALUES,&is_A_all);CHKERRQ(ierr);
308     ierr = MatGetSubMatrixUnsorted(sub_schurs->A,is_A_all,is_A_all,&A);CHKERRQ(ierr);
309     ierr = ISDestroy(&is_A_all);CHKERRQ(ierr);
310 
311     ierr = MatIsSymmetric(sub_schurs->A,0.0,&is_symmetric);CHKERRQ(ierr);
312     if (is_symmetric) {
313       ierr = MatGetFactor(A,MATSOLVERMUMPS,MAT_FACTOR_CHOLESKY,&F);CHKERRQ(ierr);
314     } else {
315       ierr = MatGetFactor(A,MATSOLVERMUMPS,MAT_FACTOR_LU,&F);CHKERRQ(ierr);
316     }
317 
318     /* subsets ordered last */
319     ierr = PetscMalloc1(local_size,&idxs_schur);CHKERRQ(ierr);
320     for (i=0;i<local_size;i++) {
321       idxs_schur[i] = n_I+i+1;
322     }
323     ierr = MatMumpsSetSchurIndices(F,local_size,idxs_schur);CHKERRQ(ierr);
324     ierr = PetscFree(idxs_schur);CHKERRQ(ierr);
325 
326     /* factorization step */
327     if (is_symmetric) {
328       ierr = MatCholeskyFactorSymbolic(F,A,NULL,NULL);CHKERRQ(ierr);
329       ierr = MatCholeskyFactorNumeric(F,A,NULL);CHKERRQ(ierr);
330     } else {
331       ierr = MatLUFactorSymbolic(F,A,NULL,NULL,NULL);CHKERRQ(ierr);
332       ierr = MatLUFactorNumeric(F,A,NULL);CHKERRQ(ierr);
333     }
334 
335     /* get explicit Schur Complement computed during numeric factorization */
336     ierr = MatMumpsGetSchurComplement(F,&S_all);CHKERRQ(ierr);
337 
338     /* free workspace */
339     ierr = MatDestroy(&A);CHKERRQ(ierr);
340     ierr = MatDestroy(&F);CHKERRQ(ierr);
341 
342 #endif
343   }
344 
345   /* TODO: just for compatibility with the previous version, needs to be fixed */
346   for (i=0;i<sub_schurs->n_subs_par;i++) {
347     PetscInt j = sub_schurs->index_parallel[i];
348     ierr = MatCreateVecs(sub_schurs->S_Ej[j],&sub_schurs->work1[j],&sub_schurs->work2[j]);CHKERRQ(ierr);
349   }
350   for (i=0;i<sub_schurs->n_subs_seq;i++) {
351      sub_schurs->work1[sub_schurs->index_sequential[i]] = 0;
352      sub_schurs->work2[sub_schurs->index_sequential[i]] = 0;
353   }
354 
355   if (!sub_schurs->n_subs_seq_g) {
356     sub_schurs->S_Ej_all = 0;
357     sub_schurs->sum_S_Ej_all = 0;
358     sub_schurs->is_Ej_all = 0;
359     PetscFunctionReturn(0);
360   }
361 
362   /* subset indices in local boundary numbering */
363   ierr = ISGlobalToLocalMappingApply(sub_schurs->BtoNmap,IS_GTOLM_DROP,local_size,all_local_idx_N+extra,&subset_size,all_local_idx_B);CHKERRQ(ierr);
364   if (subset_size != local_size) {
365     SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in sub_schurs serial (BtoNmap)! %d != %d\n",subset_size,local_size);
366   }
367   ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size,all_local_idx_B,PETSC_OWN_POINTER,&sub_schurs->is_Ej_all);CHKERRQ(ierr);
368 
369   /* Local matrix of all local Schur on subsets transposed */
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   ierr = PetscFree(nnz);CHKERRQ(ierr);
375 
376   if (!sub_schurs->use_mumps) {
377     PetscScalar *fill_vals;
378     PetscInt    *dummy_idx;
379 
380     /* Work arrays */
381     ierr = PetscMalloc1(max_subset_size,&dummy_idx);CHKERRQ(ierr);
382 
383     /* Loop on local problems end compute Schur complements explicitly */
384     local_size = 0;
385     for (i=0;i<sub_schurs->n_subs_seq;i++) {
386       Mat       S_Ej_expl;
387       PetscInt  j,lpi = sub_schurs->index_sequential[i];
388       PetscBool Sdense;
389 
390       ierr = PCBDDCComputeExplicitSchur(sub_schurs->S_Ej[lpi],&S_Ej_expl);CHKERRQ(ierr);
391       ierr = ISGetLocalSize(sub_schurs->is_AEj_B[lpi],&subset_size);CHKERRQ(ierr);
392       ierr = PetscObjectTypeCompare((PetscObject)S_Ej_expl,MATSEQDENSE,&Sdense);CHKERRQ(ierr);
393       if (Sdense) {
394         for (j=0;j<subset_size;j++) {
395           dummy_idx[j]=local_size+j;
396         }
397         ierr = MatDenseGetArray(S_Ej_expl,&fill_vals);CHKERRQ(ierr);
398         ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,fill_vals,INSERT_VALUES);CHKERRQ(ierr);
399         ierr = MatDenseRestoreArray(S_Ej_expl,&fill_vals);CHKERRQ(ierr);
400       } else {
401         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not yet implemented for sparse matrices");
402       }
403       local_size += subset_size;
404       ierr = MatDestroy(&sub_schurs->S_Ej[lpi]);CHKERRQ(ierr);
405       sub_schurs->S_Ej[lpi] = S_Ej_expl;
406     }
407     /* Stildas are not computed without mumps */
408     for (i=0;i<sub_schurs->n_subs;i++) {
409       sub_schurs->S_Ej_tilda[i] = 0;
410     }
411     ierr = PetscFree(dummy_idx);CHKERRQ(ierr);
412   } else {
413     PetscInt  *dummy_idx,n_all;
414     PetscBool is_symmetric,compute_Stilda=PETSC_TRUE;
415 
416     /* Work arrays */
417     ierr = PetscMalloc1(max_subset_size,&dummy_idx);CHKERRQ(ierr);
418 
419     /* Work arrays */
420     ierr = MatGetSize(S_all,&n_all,NULL);CHKERRQ(ierr);
421     local_size = 0;
422     ierr = MatIsSymmetric(sub_schurs->A,0.0,&is_symmetric);CHKERRQ(ierr);
423     for (i=0;i<sub_schurs->n_subs;i++) {
424       Mat S_EE;
425       IS  is_E;
426       PetscScalar *vals;
427       PetscInt j;
428 
429       ierr = ISGetLocalSize(sub_schurs->is_AEj_B[i],&subset_size);CHKERRQ(ierr);
430       ierr = ISCreateStride(PetscObjectComm((PetscObject)sub_schurs->is_I),subset_size,local_size,1,&is_E);CHKERRQ(ierr);
431       ierr = MatGetSubMatrix(S_all,is_E,is_E,MAT_INITIAL_MATRIX,&S_EE);CHKERRQ(ierr);
432       sub_schurs->S_Ej[i] = S_EE;
433       if (compute_Stilda) {
434         Mat Stilda;
435         if (sub_schurs->n_subs == 1) {
436           ierr = PetscObjectReference((PetscObject)sub_schurs->S_Ej[i]);CHKERRQ(ierr);
437           Stilda = sub_schurs->S_Ej[i];
438         } else {
439           KSP ksp;
440           PC  pc;
441           Mat S_EF,S_FE,S_FF,Stilda_impl;
442           IS  is_F;
443           PetscScalar eps=1.e-8;
444           PetscBool chop=PETSC_FALSE;
445 
446           ierr = ISComplement(is_E,0,n_all,&is_F);CHKERRQ(ierr);
447           ierr = MatGetSubMatrix(S_all,is_E,is_F,MAT_INITIAL_MATRIX,&S_EF);CHKERRQ(ierr);
448           ierr = MatGetSubMatrix(S_all,is_F,is_F,MAT_INITIAL_MATRIX,&S_FF);CHKERRQ(ierr);
449           ierr = MatGetSubMatrix(S_all,is_F,is_E,MAT_INITIAL_MATRIX,&S_FE);CHKERRQ(ierr);
450           if (chop) {
451             ierr = MatChop(S_FF,eps);CHKERRQ(ierr);
452             ierr = MatConvert(S_FF,MATAIJ,MAT_REUSE_MATRIX,&S_FF);CHKERRQ(ierr);
453           }
454           ierr = MatCreateSchurComplement(S_FF,S_FF,S_FE,S_EF,S_EE,&Stilda_impl);CHKERRQ(ierr);
455           ierr = MatSchurComplementGetKSP(Stilda_impl,&ksp);CHKERRQ(ierr);
456           ierr = KSPSetType(ksp,KSPPREONLY);CHKERRQ(ierr);
457           ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr);
458           if (is_symmetric) {
459             ierr = PCSetType(pc,PCCHOLESKY);CHKERRQ(ierr);
460           } else {
461             ierr = PCSetType(pc,PCLU);CHKERRQ(ierr);
462           }
463           if (!chop) {
464             ierr = PCFactorSetUseInPlace(pc,PETSC_TRUE);CHKERRQ(ierr);
465           } else {
466             ierr = PCFactorSetMatSolverPackage(pc,MATSOLVERMUMPS);CHKERRQ(ierr);
467           }
468           ierr = KSPSetUp(ksp);CHKERRQ(ierr);
469           ierr = PCBDDCComputeExplicitSchur(Stilda_impl,&Stilda);CHKERRQ(ierr);
470           ierr = MatDestroy(&S_FF);CHKERRQ(ierr);
471           ierr = MatDestroy(&S_FE);CHKERRQ(ierr);
472           ierr = MatDestroy(&S_EF);CHKERRQ(ierr);
473           ierr = MatDestroy(&Stilda_impl);CHKERRQ(ierr);
474           ierr = ISDestroy(&is_F);CHKERRQ(ierr);
475         }
476 /*
477         PetscViewerPushFormat(PETSC_VIEWER_STDOUT_SELF,PETSC_VIEWER_ASCII_MATLAB);
478         ierr = MatView(Stilda,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
479         PetscViewerPopFormat(PETSC_VIEWER_STDOUT_SELF);
480 */
481         sub_schurs->S_Ej_tilda[i]= Stilda;
482       } else {
483         sub_schurs->S_Ej_tilda[i] = 0;
484       }
485 
486       for (j=0;j<subset_size;j++) {
487         dummy_idx[j]=local_size+j;
488       }
489       ierr = MatDenseGetArray(S_EE,&vals);CHKERRQ(ierr);
490       ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,vals,INSERT_VALUES);CHKERRQ(ierr);
491       ierr = MatDenseRestoreArray(S_EE,&vals);CHKERRQ(ierr);
492       ierr = ISDestroy(&is_E);CHKERRQ(ierr);
493       local_size += subset_size;
494     }
495     ierr = PetscFree(dummy_idx);CHKERRQ(ierr);
496   }
497   ierr = MatAssemblyBegin(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
498   ierr = MatAssemblyEnd(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
499   /* Global matrix of all assembled Schur on subsets */
500   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);
501   ierr = ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)sub_schurs->l2gmap),1,local_size,all_local_idx_G,PETSC_COPY_VALUES,&l2gmap_subsets);CHKERRQ(ierr);
502   ierr = MatCreateIS(PetscObjectComm((PetscObject)sub_schurs->l2gmap),1,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size,l2gmap_subsets,&work_mat);CHKERRQ(ierr);
503   ierr = ISLocalToGlobalMappingDestroy(&l2gmap_subsets);CHKERRQ(ierr);
504   ierr = MatISSetLocalMat(work_mat,sub_schurs->S_Ej_all);CHKERRQ(ierr);
505   ierr = MatISGetMPIXAIJ(work_mat,MAT_INITIAL_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
506 
507   /* Get local part of (\sum_j S_Ej) */
508   ierr = PetscFree(all_local_idx_N);CHKERRQ(ierr);
509   ierr = ISCreateGeneral(PetscObjectComm((PetscObject)sub_schurs->l2gmap),local_size,all_local_idx_G,PETSC_OWN_POINTER,&temp_is);CHKERRQ(ierr);
510   ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
511   ierr = MatDestroy(&global_schur_subsets);CHKERRQ(ierr);
512 
513   /* Computation of Stildas */
514   if (S_all) {
515     PetscInt    n_all;
516     PetscScalar *vals;
517 
518     /* Work arrays */
519     ierr = MatGetSize(S_all,&n_all,NULL);CHKERRQ(ierr);
520     ierr = MatDenseGetArray(S_all,&vals);CHKERRQ(ierr);
521     ierr = MatDenseRestoreArray(S_all,&vals);CHKERRQ(ierr);
522   }
523 
524   ierr = MatDestroy(&S_all);CHKERRQ(ierr);
525   ierr = MatDestroy(&work_mat);CHKERRQ(ierr);
526   ierr = ISDestroy(&temp_is);CHKERRQ(ierr);
527   PetscFunctionReturn(0);
528 }
529 
530 #undef __FUNCT__
531 #define __FUNCT__ "PCBDDCSubSchursInit"
532 PetscErrorCode PCBDDCSubSchursInit(PCBDDCSubSchurs sub_schurs, Mat A, Mat S, IS is_I, IS is_B, PCBDDCGraph graph, ISLocalToGlobalMapping BtoNmap, PetscInt seqthreshold)
533 {
534   IS                  *faces,*edges,*all_cc;
535   PetscInt            *index_sequential,*index_parallel;
536   PetscInt            *auxlocal_sequential,*auxlocal_parallel;
537   PetscInt            *auxglobal_sequential,*auxglobal_parallel;
538   PetscInt            *auxmapping;//,*idxs;
539   PetscInt            i,max_subset_size;
540   PetscInt            n_sequential_problems,n_local_sequential_problems,n_parallel_problems,n_local_parallel_problems;
541   PetscInt            n_faces,n_edges,n_all_cc;
542   PetscBool           is_sorted;
543   PetscErrorCode  ierr;
544 
545   PetscFunctionBegin;
546   ierr = ISSorted(is_I,&is_sorted);CHKERRQ(ierr);
547   if (!is_sorted) {
548     SETERRQ(PetscObjectComm((PetscObject)is_I),PETSC_ERR_PLIB,"IS for I dofs should be shorted");
549   }
550   ierr = ISSorted(is_B,&is_sorted);CHKERRQ(ierr);
551   if (!is_sorted) {
552     SETERRQ(PetscObjectComm((PetscObject)is_B),PETSC_ERR_PLIB,"IS for B dofs should be shorted");
553   }
554 
555   /* reset any previous data */
556   ierr = PCBDDCSubSchursReset(sub_schurs);CHKERRQ(ierr);
557 
558   /* get index sets for faces and edges */
559   ierr = PCBDDCGraphGetCandidatesIS(graph,&n_faces,&faces,&n_edges,&edges,NULL);CHKERRQ(ierr);
560   n_all_cc = n_faces+n_edges;
561   ierr = PetscMalloc1(n_all_cc,&all_cc);CHKERRQ(ierr);
562   for (i=0;i<n_faces;i++) {
563     all_cc[i] = faces[i];
564   }
565   for (i=0;i<n_edges;i++) {
566     all_cc[n_faces+i] = edges[i];
567   }
568   ierr = PetscFree(faces);CHKERRQ(ierr);
569   ierr = PetscFree(edges);CHKERRQ(ierr);
570 
571   /* map interface's subsets */
572   max_subset_size = 0;
573   for (i=0;i<n_all_cc;i++) {
574     PetscInt subset_size;
575     ierr = ISGetLocalSize(all_cc[i],&subset_size);CHKERRQ(ierr);
576     max_subset_size = PetscMax(max_subset_size,subset_size);
577   }
578   ierr = PetscMalloc1(max_subset_size,&auxmapping);CHKERRQ(ierr);
579   ierr = PetscMalloc2(graph->ncc,&auxlocal_sequential,
580                       graph->ncc,&auxlocal_parallel);CHKERRQ(ierr);
581   ierr = PetscMalloc2(graph->ncc,&index_sequential,
582                       graph->ncc,&index_parallel);CHKERRQ(ierr);
583 
584   /* if threshold is negative uses all sequential problems (possibly using MUMPS) */
585   sub_schurs->use_mumps = PETSC_FALSE;
586   if (seqthreshold < 0) {
587     seqthreshold = max_subset_size;
588 #if defined(PETSC_HAVE_MUMPS)
589     sub_schurs->use_mumps = !!A;
590 #endif
591   }
592 
593 
594   /* determine which problem has to be solved in parallel or sequentially */
595   n_local_sequential_problems = 0;
596   n_local_parallel_problems = 0;
597   for (i=0;i<n_all_cc;i++) {
598     PetscInt       subset_size,j,min_loc = 0;
599     const PetscInt *idxs;
600 
601     ierr = ISGetLocalSize(all_cc[i],&subset_size);CHKERRQ(ierr);
602     ierr = ISGetIndices(all_cc[i],&idxs);CHKERRQ(ierr);
603     ierr = ISLocalToGlobalMappingApply(graph->l2gmap,subset_size,idxs,auxmapping);CHKERRQ(ierr);
604     for (j=1;j<subset_size;j++) {
605       if (auxmapping[j]<auxmapping[min_loc]) {
606         min_loc = j;
607       }
608     }
609     if (subset_size > seqthreshold) {
610       index_parallel[n_local_parallel_problems] = i;
611       auxlocal_parallel[n_local_parallel_problems] = idxs[min_loc];
612       n_local_parallel_problems++;
613     } else {
614       index_sequential[n_local_sequential_problems] = i;
615       auxlocal_sequential[n_local_sequential_problems] = idxs[min_loc];
616       n_local_sequential_problems++;
617     }
618     ierr = ISRestoreIndices(all_cc[i],&idxs);CHKERRQ(ierr);
619   }
620 
621   /* Number parallel problems */
622   auxglobal_parallel = 0;
623   ierr = PCBDDCSubsetNumbering(PetscObjectComm((PetscObject)graph->l2gmap),graph->l2gmap,n_local_parallel_problems,auxlocal_parallel,PETSC_NULL,&n_parallel_problems,&auxglobal_parallel);CHKERRQ(ierr);
624 
625   /* Number sequential problems */
626   auxglobal_sequential = 0;
627   ierr = PCBDDCSubsetNumbering(PetscObjectComm((PetscObject)graph->l2gmap),graph->l2gmap,n_local_sequential_problems,auxlocal_sequential,PETSC_NULL,&n_sequential_problems,&auxglobal_sequential);CHKERRQ(ierr);
628 
629   /* update info in sub_schurs */
630   if (sub_schurs->use_mumps && A) {
631     ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
632     sub_schurs->A = A;
633   }
634   ierr = PetscObjectReference((PetscObject)S);CHKERRQ(ierr);
635   sub_schurs->S = S;
636   ierr = PetscObjectReference((PetscObject)is_I);CHKERRQ(ierr);
637   sub_schurs->is_I = is_I;
638   ierr = PetscObjectReference((PetscObject)is_B);CHKERRQ(ierr);
639   sub_schurs->is_B = is_B;
640   ierr = PetscObjectReference((PetscObject)graph->l2gmap);CHKERRQ(ierr);
641   sub_schurs->l2gmap = graph->l2gmap;
642   ierr = PetscObjectReference((PetscObject)BtoNmap);CHKERRQ(ierr);
643   sub_schurs->BtoNmap = BtoNmap;
644   sub_schurs->n_subs_seq = n_local_sequential_problems;
645   sub_schurs->n_subs_par = n_local_parallel_problems;
646   sub_schurs->n_subs_seq_g = n_sequential_problems;
647   sub_schurs->n_subs_par_g = n_parallel_problems;
648   sub_schurs->n_subs = sub_schurs->n_subs_seq + sub_schurs->n_subs_par;
649   sub_schurs->is_subs = all_cc;
650   sub_schurs->index_sequential = index_sequential;
651   sub_schurs->index_parallel = index_parallel;
652   sub_schurs->auxglobal_sequential = auxglobal_sequential;
653   sub_schurs->auxglobal_parallel = auxglobal_parallel;
654 
655   /* free workspace */
656   ierr = PetscFree(auxmapping);CHKERRQ(ierr);
657   ierr = PetscFree2(auxlocal_sequential,auxlocal_parallel);CHKERRQ(ierr);
658 
659   PetscFunctionReturn(0);
660 }
661 
662 #undef __FUNCT__
663 #define __FUNCT__ "PCBDDCSubSchursCreate"
664 PetscErrorCode PCBDDCSubSchursCreate(PCBDDCSubSchurs *sub_schurs)
665 {
666   PCBDDCSubSchurs schurs_ctx;
667   PetscErrorCode  ierr;
668 
669   PetscFunctionBegin;
670   ierr = PetscNew(&schurs_ctx);CHKERRQ(ierr);
671   schurs_ctx->n_subs = 0;
672   *sub_schurs = schurs_ctx;
673   PetscFunctionReturn(0);
674 }
675 
676 #undef __FUNCT__
677 #define __FUNCT__ "PCBDDCSubSchursDestroy"
678 PetscErrorCode PCBDDCSubSchursDestroy(PCBDDCSubSchurs *sub_schurs)
679 {
680   PetscErrorCode ierr;
681 
682   PetscFunctionBegin;
683   ierr = PCBDDCSubSchursReset(*sub_schurs);CHKERRQ(ierr);
684   ierr = PetscFree(*sub_schurs);CHKERRQ(ierr);
685   PetscFunctionReturn(0);
686 }
687 
688 #undef __FUNCT__
689 #define __FUNCT__ "PCBDDCSubSchursReset"
690 PetscErrorCode PCBDDCSubSchursReset(PCBDDCSubSchurs sub_schurs)
691 {
692   PetscInt       i;
693   PetscErrorCode ierr;
694 
695   PetscFunctionBegin;
696   ierr = MatDestroy(&sub_schurs->A);CHKERRQ(ierr);
697   ierr = MatDestroy(&sub_schurs->S);CHKERRQ(ierr);
698   ierr = ISDestroy(&sub_schurs->is_I);CHKERRQ(ierr);
699   ierr = ISDestroy(&sub_schurs->is_B);CHKERRQ(ierr);
700   ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->l2gmap);CHKERRQ(ierr);
701   ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->BtoNmap);CHKERRQ(ierr);
702   ierr = MatDestroy(&sub_schurs->S_Ej_all);CHKERRQ(ierr);
703   ierr = MatDestroy(&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
704   ierr = ISDestroy(&sub_schurs->is_Ej_all);CHKERRQ(ierr);
705   for (i=0;i<sub_schurs->n_subs;i++) {
706     ierr = ISDestroy(&sub_schurs->is_subs[i]);CHKERRQ(ierr);
707     ierr = ISDestroy(&sub_schurs->is_AEj_B[i]);CHKERRQ(ierr);
708     ierr = MatDestroy(&sub_schurs->S_Ej[i]);CHKERRQ(ierr);
709     ierr = MatDestroy(&sub_schurs->S_Ej_tilda[i]);CHKERRQ(ierr);
710     ierr = VecDestroy(&sub_schurs->work1[i]);CHKERRQ(ierr);
711     ierr = VecDestroy(&sub_schurs->work2[i]);CHKERRQ(ierr);
712   }
713   ierr = ISDestroy(&sub_schurs->is_I_layer);CHKERRQ(ierr);
714   if (sub_schurs->n_subs) {
715     ierr = PetscFree(sub_schurs->is_subs);CHKERRQ(ierr);
716     ierr = PetscFree5(sub_schurs->is_AEj_B,sub_schurs->S_Ej,sub_schurs->S_Ej_tilda,sub_schurs->work1,sub_schurs->work2);CHKERRQ(ierr);
717     ierr = PetscFree2(sub_schurs->index_sequential,sub_schurs->index_parallel);CHKERRQ(ierr);
718     ierr = PetscFree(sub_schurs->auxglobal_sequential);CHKERRQ(ierr);
719     ierr = PetscFree(sub_schurs->auxglobal_parallel);CHKERRQ(ierr);
720   }
721   sub_schurs->n_subs = 0;
722   PetscFunctionReturn(0);
723 }
724 
725 #undef __FUNCT__
726 #define __FUNCT__ "PCBDDCAdjGetNextLayer_Private"
727 PETSC_STATIC_INLINE PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt* queue_tip,PetscInt n_prev,PetscBT touched,PetscInt* xadj,PetscInt* adjncy,PetscInt* n_added)
728 {
729   PetscInt       i,j,n;
730   PetscErrorCode ierr;
731 
732   PetscFunctionBegin;
733   n = 0;
734   for (i=-n_prev;i<0;i++) {
735     PetscInt start_dof = queue_tip[i];
736     for (j=xadj[start_dof];j<xadj[start_dof+1];j++) {
737       PetscInt dof = adjncy[j];
738       if (!PetscBTLookup(touched,dof)) {
739         ierr = PetscBTSet(touched,dof);CHKERRQ(ierr);
740         queue_tip[n] = dof;
741         n++;
742       }
743     }
744   }
745   *n_added = n;
746   PetscFunctionReturn(0);
747 }
748