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