xref: /petsc/src/ksp/pc/impls/bddc/bddcschurs.c (revision d648f8588b208cff93aa1a93856aa59d269ed7cc)
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_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   extra = 0;
216   if (sub_schurs->use_mumps) {
217     ierr = ISGetLocalSize(sub_schurs->is_I_layer,&extra);CHKERRQ(ierr);
218   }
219   ierr = PetscMalloc1(local_size+extra,&all_local_idx_N);CHKERRQ(ierr);
220   if (extra) {
221     const PetscInt *idxs;
222     ierr = ISGetIndices(sub_schurs->is_I_layer,&idxs);CHKERRQ(ierr);
223     ierr = PetscMemcpy(all_local_idx_N,idxs,extra*sizeof(PetscInt));CHKERRQ(ierr);
224     ierr = ISRestoreIndices(sub_schurs->is_I_layer,&idxs);CHKERRQ(ierr);
225   }
226   ierr = PetscMalloc1(local_size,&nnz);CHKERRQ(ierr);
227 
228   /* Get local indices in local numbering */
229   local_size = 0;
230   for (i=0;i<sub_schurs->n_subs_seq;i++) {
231     PetscInt j;
232     const    PetscInt *idxs;
233 
234     PetscInt local_problem_index = sub_schurs->index_sequential[i];
235     ierr = ISGetLocalSize(sub_schurs->is_subs[local_problem_index],&subset_size);CHKERRQ(ierr);
236     ierr = ISGetIndices(sub_schurs->is_subs[local_problem_index],&idxs);CHKERRQ(ierr);
237     /* subset indices in local numbering */
238     ierr = PetscMemcpy(all_local_idx_N+local_size+extra,idxs,subset_size*sizeof(PetscInt));CHKERRQ(ierr);
239     ierr = ISRestoreIndices(sub_schurs->is_subs[local_problem_index],&idxs);CHKERRQ(ierr);
240     for (j=0;j<subset_size;j++) nnz[local_size+j] = subset_size;
241     local_size += subset_size;
242   }
243 
244   S_all = NULL;
245   sub_schurs->is_hermitian = PETSC_FALSE;
246   sub_schurs->is_posdef = PETSC_FALSE;
247   if (sub_schurs->A) {
248     ierr = MatIsHermitian(sub_schurs->A,0.0,&sub_schurs->is_hermitian);CHKERRQ(ierr);
249     if (sub_schurs->is_hermitian) {
250       PetscScalar val;
251       Vec         vec1,vec2;
252 
253       ierr = MatCreateVecs(sub_schurs->A,&vec1,&vec2);CHKERRQ(ierr);
254       ierr = VecSetRandom(vec1,NULL);
255       ierr = VecCopy(vec1,vec2);CHKERRQ(ierr);
256       ierr = MatMult(sub_schurs->A,vec2,vec1);CHKERRQ(ierr);
257       ierr = VecDot(vec1,vec2,&val);CHKERRQ(ierr);
258       if (PetscRealPart(val) > 0. && PetscImaginaryPart(val) == 0.) sub_schurs->is_posdef = PETSC_TRUE;
259       ierr = VecDestroy(&vec1);CHKERRQ(ierr);
260       ierr = VecDestroy(&vec2);CHKERRQ(ierr);
261     }
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,&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     PetscInt *all_local_idx_B;
391 
392     ierr = PetscMalloc1(local_size,&all_local_idx_B);CHKERRQ(ierr);
393     ierr = ISGlobalToLocalMappingApply(sub_schurs->BtoNmap,IS_GTOLM_DROP,local_size,all_local_idx_N+extra,&subset_size,all_local_idx_B);CHKERRQ(ierr);
394     if (subset_size != local_size) {
395       SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in sub_schurs serial (BtoNmap)! %d != %d\n",subset_size,local_size);
396     }
397     ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size,all_local_idx_B,PETSC_OWN_POINTER,&sub_schurs->is_Ej_all);CHKERRQ(ierr);
398   }
399 
400   /* Local matrix of all local Schur on subsets transposed */
401   if (!sub_schurs->S_Ej_all) {
402     ierr = MatCreate(PETSC_COMM_SELF,&sub_schurs->S_Ej_all);CHKERRQ(ierr);
403     ierr = MatSetSizes(sub_schurs->S_Ej_all,PETSC_DECIDE,PETSC_DECIDE,local_size,local_size);CHKERRQ(ierr);
404     ierr = MatSetType(sub_schurs->S_Ej_all,MATAIJ);CHKERRQ(ierr);
405     ierr = MatSeqAIJSetPreallocation(sub_schurs->S_Ej_all,0,nnz);CHKERRQ(ierr);
406   } else {
407     ierr = MatZeroEntries(sub_schurs->S_Ej_all);CHKERRQ(ierr);
408   }
409 
410   S_Ej_tilda_all = 0;
411   S_Ej_inv_all = 0;
412   Bwork = NULL;
413   pivots = NULL;
414   if (sub_schurs->n_subs && deluxe && compute_Stilda && !sub_schurs->is_hermitian) { /* workspace needed only for GETRI */
415     PetscScalar lwork;
416 
417     B_lwork = -1;
418     ierr = PetscBLASIntCast(max_subset_size_seq,&B_N);CHKERRQ(ierr);
419     ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
420     PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,Bwork,&B_N,pivots,&lwork,&B_lwork,&B_ierr));
421     ierr = PetscFPTrapPop();CHKERRQ(ierr);
422     if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to GETRI Lapack routine %d",(int)B_ierr);
423     ierr = PetscBLASIntCast((PetscInt)lwork,&B_lwork);CHKERRQ(ierr);
424     ierr = PetscMalloc2(B_lwork,&Bwork,max_subset_size_seq,&pivots);CHKERRQ(ierr);
425   }
426 
427   ierr = PetscBTMemzero(sub_schurs->n_subs,sub_schurs->computed_Stilda_subs);CHKERRQ(ierr);
428   if (!sub_schurs->use_mumps) {
429     if (sub_schurs->n_subs_seq) {
430       PetscScalar *work;
431       PetscInt    *dummy_idx;
432 
433       /* Work arrays */
434       ierr = PetscMalloc2(max_subset_size_seq,&dummy_idx,max_subset_size_seq*max_subset_size_seq,&work);CHKERRQ(ierr);
435 
436       /* Loop on local problems end compute Schur complements explicitly */
437       local_size = 0;
438       for (i=0;i<sub_schurs->n_subs_seq;i++) {
439         Mat       S_Ej_expl;
440         PetscInt  j,lpi = sub_schurs->index_sequential[i];
441         PetscBool Sdense;
442 
443         ierr = ISGetLocalSize(sub_schurs->is_subs[lpi],&subset_size);CHKERRQ(ierr);
444         ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&S_Ej_expl);CHKERRQ(ierr);
445         ierr = PCBDDCComputeExplicitSchur(sub_schurs->S_Ej[lpi],sub_schurs->is_hermitian,MAT_REUSE_MATRIX,&S_Ej_expl);CHKERRQ(ierr);
446         ierr = PetscObjectTypeCompare((PetscObject)S_Ej_expl,MATSEQDENSE,&Sdense);CHKERRQ(ierr);
447         if (Sdense) {
448           for (j=0;j<subset_size;j++) {
449             dummy_idx[j]=local_size+j;
450           }
451           ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
452         } else {
453           SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not yet implemented for sparse matrices");
454         }
455         local_size += subset_size;
456         ierr = MatDestroy(&sub_schurs->S_Ej[lpi]);CHKERRQ(ierr);
457         sub_schurs->S_Ej[lpi] = S_Ej_expl;
458       }
459       ierr = PetscFree2(dummy_idx,work);CHKERRQ(ierr);
460     }
461   } else {
462     PetscInt    *dummy_idx,n_all;
463     PetscScalar *work;
464 
465     if (compute_Stilda) {
466       ierr = MatCreate(PETSC_COMM_SELF,&S_Ej_tilda_all);CHKERRQ(ierr);
467       ierr = MatSetSizes(S_Ej_tilda_all,PETSC_DECIDE,PETSC_DECIDE,local_size,local_size);CHKERRQ(ierr);
468       ierr = MatSetType(S_Ej_tilda_all,MATAIJ);CHKERRQ(ierr);
469       ierr = MatSeqAIJSetPreallocation(S_Ej_tilda_all,0,nnz);CHKERRQ(ierr);
470       if (deluxe) {
471         ierr = MatCreate(PETSC_COMM_SELF,&S_Ej_inv_all);CHKERRQ(ierr);
472         ierr = MatSetSizes(S_Ej_inv_all,PETSC_DECIDE,PETSC_DECIDE,local_size,local_size);CHKERRQ(ierr);
473         ierr = MatSetType(S_Ej_inv_all,MATAIJ);CHKERRQ(ierr);
474         ierr = MatSeqAIJSetPreallocation(S_Ej_inv_all,0,nnz);CHKERRQ(ierr);
475       }
476     }
477 
478     if (S_all) { /* multilevel */
479       ierr = MatGetSize(S_all,&n_all,NULL);CHKERRQ(ierr);
480     } else {
481       n_all = 0;
482     }
483     local_size = 0;
484 
485     /* Work arrays */
486     if (compute_Stilda) {
487       ierr = PetscMalloc2(max_subset_size_seq,&dummy_idx,n_all*n_all,&work);CHKERRQ(ierr);
488     } else {
489       ierr = PetscMalloc2(max_subset_size_seq,&dummy_idx,0,&work);CHKERRQ(ierr);
490     }
491 
492     for (i=0;i<sub_schurs->n_subs;i++) {
493       IS  is_E;
494       PetscScalar *vals;
495       PetscInt j;
496 
497       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
498       for (j=0;j<subset_size;j++) {
499         dummy_idx[j]=local_size+j;
500       }
501       ierr = ISCreateStride(PetscObjectComm((PetscObject)sub_schurs->is_I),subset_size,local_size,1,&is_E);CHKERRQ(ierr);
502       if (sub_schurs->S_Ej[i]) {
503         ierr = MatGetSubMatrix(S_all,is_E,is_E,MAT_REUSE_MATRIX,&sub_schurs->S_Ej[i]);CHKERRQ(ierr);
504       } else {
505         ierr = MatGetSubMatrix(S_all,is_E,is_E,MAT_INITIAL_MATRIX,&sub_schurs->S_Ej[i]);CHKERRQ(ierr);
506       }
507       ierr = MatDenseGetArray(sub_schurs->S_Ej[i],&vals);CHKERRQ(ierr);
508       ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,vals,INSERT_VALUES);CHKERRQ(ierr);
509       ierr = MatDenseRestoreArray(sub_schurs->S_Ej[i],&vals);CHKERRQ(ierr);
510 
511       if (compute_Stilda && ((PetscBTLookup(sub_schurs->is_edge,i) && use_edges) || use_faces)) {
512         Mat Stilda;
513 
514         ierr = PetscBTSet(sub_schurs->computed_Stilda_subs,i);CHKERRQ(ierr);
515         if (sub_schurs->n_subs == 1) {
516           ierr = PetscObjectReference((PetscObject)sub_schurs->S_Ej[i]);CHKERRQ(ierr);
517           Stilda = sub_schurs->S_Ej[i];
518         } else {
519           KSP         ksp;
520           PC          pc;
521           Mat         S_EF,S_FE,S_FF,Stilda_impl;
522           IS          is_F;
523           PetscScalar eps=1.e-8;
524           PetscInt    nc,cum;
525           PetscBool   chop=PETSC_FALSE;
526 
527           ierr = ISComplement(is_E,0,n_all,&is_F);CHKERRQ(ierr);
528           nc   = n_all - subset_size;
529           cum  = 0;
530           ierr = MatCreateSeqDense(PETSC_COMM_SELF,nc,nc,work+cum,&S_FF);CHKERRQ(ierr);
531           cum += nc*nc;
532           ierr = MatCreateSeqDense(PETSC_COMM_SELF,nc,subset_size,work+cum,&S_FE);CHKERRQ(ierr);
533           cum += nc*subset_size;
534           ierr = MatGetSubMatrix(S_all,is_F,is_F,MAT_REUSE_MATRIX,&S_FF);CHKERRQ(ierr);
535           ierr = MatGetSubMatrix(S_all,is_F,is_E,MAT_REUSE_MATRIX,&S_FE);CHKERRQ(ierr);
536           if (sub_schurs->is_hermitian) { /* just need a placeholder for S_EF */
537             ierr = MatCreateTranspose(S_FE,&S_EF);CHKERRQ(ierr);
538           } else {
539             ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,nc,work+cum,&S_EF);CHKERRQ(ierr);
540             cum += nc*subset_size;
541             ierr = MatGetSubMatrix(S_all,is_E,is_F,MAT_REUSE_MATRIX,&S_EF);CHKERRQ(ierr);
542           }
543           ierr = ISDestroy(&is_F);CHKERRQ(ierr);
544           if (chop) {
545             ierr = MatChop(S_FF,eps);CHKERRQ(ierr);
546             ierr = MatConvert(S_FF,MATAIJ,MAT_REUSE_MATRIX,&S_FF);CHKERRQ(ierr);
547           }
548           ierr = MatCreateSchurComplement(S_FF,S_FF,S_FE,S_EF,sub_schurs->S_Ej[i],&Stilda_impl);CHKERRQ(ierr);
549           ierr = MatDestroy(&S_FF);CHKERRQ(ierr);
550           ierr = MatDestroy(&S_FE);CHKERRQ(ierr);
551           ierr = MatDestroy(&S_EF);CHKERRQ(ierr);
552           ierr = MatSchurComplementGetKSP(Stilda_impl,&ksp);CHKERRQ(ierr);
553           ierr = KSPSetType(ksp,KSPPREONLY);CHKERRQ(ierr);
554           ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr);
555           if (sub_schurs->is_hermitian && sub_schurs->is_posdef) {
556             ierr = PCSetType(pc,PCCHOLESKY);CHKERRQ(ierr);
557           } else {
558             ierr = PCSetType(pc,PCLU);CHKERRQ(ierr);
559           }
560           if (!chop) {
561             ierr = PCFactorSetUseInPlace(pc,PETSC_TRUE);CHKERRQ(ierr);
562           } else {
563             ierr = PCFactorSetMatSolverPackage(pc,MATSOLVERMUMPS);CHKERRQ(ierr);
564           }
565           ierr = KSPSetUp(ksp);CHKERRQ(ierr);
566           ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work+cum,&Stilda);CHKERRQ(ierr);
567           ierr = PCBDDCComputeExplicitSchur(Stilda_impl,sub_schurs->is_hermitian,MAT_REUSE_MATRIX,&Stilda);CHKERRQ(ierr);
568           ierr = MatDestroy(&Stilda_impl);CHKERRQ(ierr);
569         }
570 /*
571         PetscViewerPushFormat(PETSC_VIEWER_STDOUT_SELF,PETSC_VIEWER_ASCII_MATLAB);
572         ierr = MatView(Stilda,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
573         PetscViewerPopFormat(PETSC_VIEWER_STDOUT_SELF);
574 */
575         ierr = MatDenseGetArray(Stilda,&vals);CHKERRQ(ierr);
576         if (deluxe) { /* when using deluxe scaling, we need (S_1^-1+S_2^-1)^-1 */
577           PetscScalar *vals2;
578 
579           ierr = MatDenseGetArray(sub_schurs->S_Ej[i],&vals2);CHKERRQ(ierr);
580           ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr);
581           ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
582           if (invert_Stildas) { /* Stildas can be singular */
583             if (!sub_schurs->is_hermitian) {
584               PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,vals,&B_N,pivots,&B_ierr));
585               if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
586               PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,vals,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
587               if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
588             } else {
589               PetscInt j,k;
590 
591               PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,vals,&B_N,&B_ierr));
592               if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
593               PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,vals,&B_N,&B_ierr));
594               if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
595               for (j=0;j<B_N;j++) {
596                 for (k=j+1;k<B_N;k++) {
597                   vals[k*B_N+j] = vals[j*B_N+k];
598                 }
599               }
600             }
601           }
602           if (!sub_schurs->is_hermitian) {
603             PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,vals2,&B_N,pivots,&B_ierr));
604             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
605             PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,vals2,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
606             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
607           } else {
608             PetscInt j,k;
609 
610             PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,vals2,&B_N,&B_ierr));
611             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
612             PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,vals2,&B_N,&B_ierr));
613             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
614             for (j=0;j<B_N;j++) {
615               for (k=j+1;k<B_N;k++) {
616                 vals2[k*B_N+j] = vals2[j*B_N+k];
617               }
618             }
619           }
620           ierr = PetscFPTrapPop();CHKERRQ(ierr);
621           ierr = MatSetValues(S_Ej_inv_all,subset_size,dummy_idx,subset_size,dummy_idx,vals2,INSERT_VALUES);CHKERRQ(ierr);
622           ierr = MatDenseRestoreArray(sub_schurs->S_Ej[i],&vals2);CHKERRQ(ierr);
623         }
624         ierr = MatSetValues(S_Ej_tilda_all,subset_size,dummy_idx,subset_size,dummy_idx,vals,INSERT_VALUES);CHKERRQ(ierr);
625         ierr = MatDenseRestoreArray(Stilda,&vals);CHKERRQ(ierr);
626         ierr = MatDestroy(&Stilda);CHKERRQ(ierr);
627       }
628       ierr = ISDestroy(&is_E);CHKERRQ(ierr);
629       local_size += subset_size;
630     }
631     ierr = PetscFree2(dummy_idx,work);CHKERRQ(ierr);
632   }
633   ierr = PetscFree(nnz);CHKERRQ(ierr);
634   ierr = MatDestroy(&S_all);CHKERRQ(ierr);
635   ierr = MatAssemblyBegin(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
636   ierr = MatAssemblyEnd(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
637   if (S_Ej_tilda_all) {
638     ierr = MatAssemblyBegin(S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
639     ierr = MatAssemblyEnd(S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
640   }
641   if (S_Ej_inv_all) {
642     ierr = MatAssemblyBegin(S_Ej_inv_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
643     ierr = MatAssemblyEnd(S_Ej_inv_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
644   }
645 
646   /* Global matrix of all assembled Schur on subsets */
647   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);
648   ierr = ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)sub_schurs->l2gmap),1,local_size,all_local_idx_G,PETSC_COPY_VALUES,&l2gmap_subsets);CHKERRQ(ierr);
649   ierr = MatCreateIS(PetscObjectComm((PetscObject)sub_schurs->l2gmap),1,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size,l2gmap_subsets,&work_mat);CHKERRQ(ierr);
650   ierr = ISLocalToGlobalMappingDestroy(&l2gmap_subsets);CHKERRQ(ierr);
651   ierr = MatISSetLocalMat(work_mat,sub_schurs->S_Ej_all);CHKERRQ(ierr);
652   ierr = MatCreate(PetscObjectComm((PetscObject)work_mat),&global_schur_subsets);CHKERRQ(ierr);
653   ierr = MatSetSizes(global_schur_subsets,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size);CHKERRQ(ierr);
654   ierr = MatSetType(global_schur_subsets,MATMPIAIJ);CHKERRQ(ierr);
655   ierr = MatISSetMPIXAIJPreallocation_Private(work_mat,global_schur_subsets,PETSC_TRUE);CHKERRQ(ierr);
656   ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
657   /* Get local part of (\sum_j S_Ej) */
658   ierr = PetscFree(all_local_idx_N);CHKERRQ(ierr);
659   ierr = ISCreateGeneral(PetscObjectComm((PetscObject)sub_schurs->l2gmap),local_size,all_local_idx_G,PETSC_OWN_POINTER,&temp_is);CHKERRQ(ierr);
660   ierr = MatDestroy(&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
661   ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
662 
663   if (S_Ej_tilda_all) {
664     ierr = MatISSetLocalMat(work_mat,S_Ej_tilda_all);CHKERRQ(ierr);
665     ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
666     ierr = MatDestroy(&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
667     ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
668   }
669   if (S_Ej_inv_all) {
670     PetscInt    cum;
671     PetscScalar *array,*array2;
672     ierr = MatISSetLocalMat(work_mat,S_Ej_inv_all);CHKERRQ(ierr);
673     ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
674     ierr = MatDestroy(&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
675     ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
676     /* invert blocks of sum_S_Ej_inv_all */
677     ierr = MatSeqAIJGetArray(sub_schurs->sum_S_Ej_inv_all,&array);CHKERRQ(ierr);
678     ierr = MatSeqAIJGetArray(sub_schurs->sum_S_Ej_tilda_all,&array2);CHKERRQ(ierr);
679     cum = 0;
680     for (i=0;i<sub_schurs->n_subs;i++) {
681       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
682       if (PetscBTLookup(sub_schurs->computed_Stilda_subs,i)) {
683         ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr);
684         ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
685         if (!sub_schurs->is_hermitian) {
686           PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,array+cum,&B_N,pivots,&B_ierr));
687           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
688           PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,array+cum,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
689           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
690         } else {
691           PetscInt j,k;
692 
693           PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,array+cum,&B_N,&B_ierr));
694           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
695           PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,array+cum,&B_N,&B_ierr));
696           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
697           for (j=0;j<B_N;j++) {
698             for (k=j+1;k<B_N;k++) {
699               array[k*B_N+j+cum] = array[j*B_N+k+cum];
700             }
701           }
702         }
703         if (invert_Stildas) { /* Stildas can be singular */
704           if (!sub_schurs->is_hermitian) {
705             PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,array2+cum,&B_N,pivots,&B_ierr));
706             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
707             PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,array2+cum,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
708             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
709           } else {
710             PetscInt j,k;
711 
712             PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,array2+cum,&B_N,&B_ierr));
713             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
714             PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,array2+cum,&B_N,&B_ierr));
715             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
716             for (j=0;j<B_N;j++) {
717               for (k=j+1;k<B_N;k++) {
718                 array2[k*B_N+j+cum] = array2[j*B_N+k+cum];
719               }
720             }
721           }
722         }
723         ierr = PetscFPTrapPop();CHKERRQ(ierr);
724       }
725       cum += subset_size*subset_size;
726     }
727     ierr = MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_inv_all,&array);CHKERRQ(ierr);
728     ierr = MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_tilda_all,&array2);CHKERRQ(ierr);
729   }
730 
731   ierr = PetscFree2(Bwork,pivots);CHKERRQ(ierr);
732   ierr = MatDestroy(&global_schur_subsets);CHKERRQ(ierr);
733   ierr = MatDestroy(&S_Ej_tilda_all);CHKERRQ(ierr);
734   ierr = MatDestroy(&S_Ej_inv_all);CHKERRQ(ierr);
735   ierr = MatDestroy(&work_mat);CHKERRQ(ierr);
736   ierr = ISDestroy(&temp_is);CHKERRQ(ierr);
737   PetscFunctionReturn(0);
738 }
739 
740 #undef __FUNCT__
741 #define __FUNCT__ "PCBDDCSubSchursInit"
742 PetscErrorCode PCBDDCSubSchursInit(PCBDDCSubSchurs sub_schurs, Mat A, Mat S, IS is_I, IS is_B, PCBDDCGraph graph, ISLocalToGlobalMapping BtoNmap, PetscInt seqthreshold)
743 {
744   IS                  *faces,*edges,*all_cc,vertices;
745   PetscInt            *index_sequential,*index_parallel;
746   PetscInt            *auxlocal_sequential,*auxlocal_parallel;
747   PetscInt            *auxglobal_sequential,*auxglobal_parallel;
748   PetscInt            *auxmapping;
749   PetscInt            i,max_subset_size;
750   PetscInt            n_sequential_problems,n_local_sequential_problems,n_parallel_problems,n_local_parallel_problems;
751   PetscInt            n_faces,n_edges,n_all_cc;
752   PetscBool           is_sorted;
753   PetscErrorCode  ierr;
754 
755   PetscFunctionBegin;
756   ierr = ISSorted(is_I,&is_sorted);CHKERRQ(ierr);
757   if (!is_sorted) {
758     SETERRQ(PetscObjectComm((PetscObject)is_I),PETSC_ERR_PLIB,"IS for I dofs should be shorted");
759   }
760   ierr = ISSorted(is_B,&is_sorted);CHKERRQ(ierr);
761   if (!is_sorted) {
762     SETERRQ(PetscObjectComm((PetscObject)is_B),PETSC_ERR_PLIB,"IS for B dofs should be shorted");
763   }
764 
765   /* reset any previous data */
766   ierr = PCBDDCSubSchursReset(sub_schurs);CHKERRQ(ierr);
767 
768   /* get index sets for faces and edges */
769   ierr = PCBDDCGraphGetCandidatesIS(graph,&n_faces,&faces,&n_edges,&edges,&vertices);CHKERRQ(ierr);
770   n_all_cc = n_faces+n_edges;
771   ierr = PetscBTCreate(n_all_cc,&sub_schurs->is_edge);CHKERRQ(ierr);
772   ierr = PetscBTCreate(n_all_cc,&sub_schurs->computed_Stilda_subs);CHKERRQ(ierr);
773   ierr = PetscMalloc1(n_all_cc,&all_cc);CHKERRQ(ierr);
774   for (i=0;i<n_faces;i++) {
775     all_cc[i] = faces[i];
776   }
777   for (i=0;i<n_edges;i++) {
778     all_cc[n_faces+i] = edges[i];
779     ierr = PetscBTSet(sub_schurs->is_edge,n_faces+i);CHKERRQ(ierr);
780   }
781   ierr = PetscFree(faces);CHKERRQ(ierr);
782   ierr = PetscFree(edges);CHKERRQ(ierr);
783 
784   /* map interface's subsets */
785   max_subset_size = 0;
786   for (i=0;i<n_all_cc;i++) {
787     PetscInt subset_size;
788     ierr = ISGetLocalSize(all_cc[i],&subset_size);CHKERRQ(ierr);
789     max_subset_size = PetscMax(max_subset_size,subset_size);
790   }
791   ierr = PetscMalloc1(max_subset_size,&auxmapping);CHKERRQ(ierr);
792   ierr = PetscMalloc2(graph->ncc,&auxlocal_sequential,
793                       graph->ncc,&auxlocal_parallel);CHKERRQ(ierr);
794   ierr = PetscMalloc2(graph->ncc,&index_sequential,
795                       graph->ncc,&index_parallel);CHKERRQ(ierr);
796 
797   /* if threshold is negative uses all sequential problems (possibly using MUMPS) */
798   sub_schurs->use_mumps = PETSC_FALSE;
799   if (seqthreshold < 0) {
800     seqthreshold = max_subset_size;
801 #if defined(PETSC_HAVE_MUMPS)
802     sub_schurs->use_mumps = !!A;
803 #endif
804   }
805 
806 
807   /* determine which problem has to be solved in parallel or sequentially */
808   n_local_sequential_problems = 0;
809   n_local_parallel_problems = 0;
810   for (i=0;i<n_all_cc;i++) {
811     PetscInt       subset_size,j,min_loc = 0;
812     const PetscInt *idxs;
813 
814     ierr = ISGetLocalSize(all_cc[i],&subset_size);CHKERRQ(ierr);
815     ierr = ISGetIndices(all_cc[i],&idxs);CHKERRQ(ierr);
816     ierr = ISLocalToGlobalMappingApply(graph->l2gmap,subset_size,idxs,auxmapping);CHKERRQ(ierr);
817     for (j=1;j<subset_size;j++) {
818       if (auxmapping[j]<auxmapping[min_loc]) {
819         min_loc = j;
820       }
821     }
822     if (subset_size > seqthreshold) {
823       index_parallel[n_local_parallel_problems] = i;
824       auxlocal_parallel[n_local_parallel_problems] = idxs[min_loc];
825       n_local_parallel_problems++;
826     } else {
827       index_sequential[n_local_sequential_problems] = i;
828       auxlocal_sequential[n_local_sequential_problems] = idxs[min_loc];
829       n_local_sequential_problems++;
830     }
831     ierr = ISRestoreIndices(all_cc[i],&idxs);CHKERRQ(ierr);
832   }
833 
834   /* Number parallel problems */
835   auxglobal_parallel = 0;
836   ierr = PCBDDCSubsetNumbering(PetscObjectComm((PetscObject)graph->l2gmap),graph->l2gmap,n_local_parallel_problems,auxlocal_parallel,PETSC_NULL,&n_parallel_problems,&auxglobal_parallel);CHKERRQ(ierr);
837 
838   /* Number sequential problems */
839   auxglobal_sequential = 0;
840   ierr = PCBDDCSubsetNumbering(PetscObjectComm((PetscObject)graph->l2gmap),graph->l2gmap,n_local_sequential_problems,auxlocal_sequential,PETSC_NULL,&n_sequential_problems,&auxglobal_sequential);CHKERRQ(ierr);
841 
842   /* update info in sub_schurs */
843   if (A) {
844     PetscBool isseqaij;
845 
846     ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
847     if (isseqaij) {
848       ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
849       sub_schurs->A = A;
850     } else { /* SeqBAIJ matrices does not support symmetry checking, SeqSBAIJ does not support MatPermute */
851       ierr = MatConvert(A,MATSEQAIJ,MAT_INITIAL_MATRIX,&sub_schurs->A);CHKERRQ(ierr);
852     }
853   }
854   ierr = PetscObjectReference((PetscObject)S);CHKERRQ(ierr);
855   sub_schurs->S = S;
856   ierr = PetscObjectReference((PetscObject)is_I);CHKERRQ(ierr);
857   sub_schurs->is_I = is_I;
858   ierr = PetscObjectReference((PetscObject)is_B);CHKERRQ(ierr);
859   sub_schurs->is_B = is_B;
860   ierr = PetscObjectReference((PetscObject)graph->l2gmap);CHKERRQ(ierr);
861   sub_schurs->l2gmap = graph->l2gmap;
862   ierr = PetscObjectReference((PetscObject)BtoNmap);CHKERRQ(ierr);
863   sub_schurs->BtoNmap = BtoNmap;
864   sub_schurs->n_subs_seq = n_local_sequential_problems;
865   sub_schurs->n_subs_par = n_local_parallel_problems;
866   sub_schurs->n_subs_seq_g = n_sequential_problems;
867   sub_schurs->n_subs_par_g = n_parallel_problems;
868   sub_schurs->n_subs = sub_schurs->n_subs_seq + sub_schurs->n_subs_par;
869   sub_schurs->is_subs = all_cc;
870   if (!sub_schurs->use_mumps) { /* for adaptive selection */
871     for (i=0;i<sub_schurs->n_subs;i++) {
872       ierr = ISSort(sub_schurs->is_subs[i]);CHKERRQ(ierr);
873     }
874   }
875   sub_schurs->is_Ej_com = vertices;
876   sub_schurs->index_sequential = index_sequential;
877   sub_schurs->index_parallel = index_parallel;
878   sub_schurs->auxglobal_sequential = auxglobal_sequential;
879   sub_schurs->auxglobal_parallel = auxglobal_parallel;
880 
881 
882   /* allocate space for schur complements */
883   ierr = PetscCalloc1(sub_schurs->n_subs,&sub_schurs->S_Ej);CHKERRQ(ierr);
884   sub_schurs->S_Ej_all = NULL;
885   sub_schurs->sum_S_Ej_all = NULL;
886   sub_schurs->sum_S_Ej_inv_all = NULL;
887   sub_schurs->sum_S_Ej_tilda_all = NULL;
888   sub_schurs->is_Ej_all = NULL;
889 
890   /* free workspace */
891   ierr = PetscFree(auxmapping);CHKERRQ(ierr);
892   ierr = PetscFree2(auxlocal_sequential,auxlocal_parallel);CHKERRQ(ierr);
893 
894   PetscFunctionReturn(0);
895 }
896 
897 #undef __FUNCT__
898 #define __FUNCT__ "PCBDDCSubSchursCreate"
899 PetscErrorCode PCBDDCSubSchursCreate(PCBDDCSubSchurs *sub_schurs)
900 {
901   PCBDDCSubSchurs schurs_ctx;
902   PetscErrorCode  ierr;
903 
904   PetscFunctionBegin;
905   ierr = PetscNew(&schurs_ctx);CHKERRQ(ierr);
906   schurs_ctx->n_subs = 0;
907   *sub_schurs = schurs_ctx;
908   PetscFunctionReturn(0);
909 }
910 
911 #undef __FUNCT__
912 #define __FUNCT__ "PCBDDCSubSchursDestroy"
913 PetscErrorCode PCBDDCSubSchursDestroy(PCBDDCSubSchurs *sub_schurs)
914 {
915   PetscErrorCode ierr;
916 
917   PetscFunctionBegin;
918   ierr = PCBDDCSubSchursReset(*sub_schurs);CHKERRQ(ierr);
919   ierr = PetscFree(*sub_schurs);CHKERRQ(ierr);
920   PetscFunctionReturn(0);
921 }
922 
923 #undef __FUNCT__
924 #define __FUNCT__ "PCBDDCSubSchursReset"
925 PetscErrorCode PCBDDCSubSchursReset(PCBDDCSubSchurs sub_schurs)
926 {
927   PetscInt       i;
928   PetscErrorCode ierr;
929 
930   PetscFunctionBegin;
931   ierr = MatDestroy(&sub_schurs->A);CHKERRQ(ierr);
932   ierr = MatDestroy(&sub_schurs->S);CHKERRQ(ierr);
933   ierr = ISDestroy(&sub_schurs->is_I);CHKERRQ(ierr);
934   ierr = ISDestroy(&sub_schurs->is_B);CHKERRQ(ierr);
935   ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->l2gmap);CHKERRQ(ierr);
936   ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->BtoNmap);CHKERRQ(ierr);
937   ierr = MatDestroy(&sub_schurs->S_Ej_all);CHKERRQ(ierr);
938   ierr = MatDestroy(&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
939   ierr = MatDestroy(&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
940   ierr = MatDestroy(&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
941   ierr = ISDestroy(&sub_schurs->is_Ej_all);CHKERRQ(ierr);
942   ierr = ISDestroy(&sub_schurs->is_Ej_com);CHKERRQ(ierr);
943   ierr = PetscBTDestroy(&sub_schurs->is_edge);CHKERRQ(ierr);
944   ierr = PetscBTDestroy(&sub_schurs->computed_Stilda_subs);CHKERRQ(ierr);
945   for (i=0;i<sub_schurs->n_subs;i++) {
946     ierr = ISDestroy(&sub_schurs->is_subs[i]);CHKERRQ(ierr);
947     ierr = MatDestroy(&sub_schurs->S_Ej[i]);CHKERRQ(ierr);
948   }
949   ierr = ISDestroy(&sub_schurs->is_I_layer);CHKERRQ(ierr);
950   if (sub_schurs->n_subs) {
951     ierr = PetscFree(sub_schurs->is_subs);CHKERRQ(ierr);
952     ierr = PetscFree(sub_schurs->S_Ej);CHKERRQ(ierr);
953   }
954   ierr = PetscFree2(sub_schurs->index_sequential,sub_schurs->index_parallel);CHKERRQ(ierr);
955   ierr = PetscFree(sub_schurs->auxglobal_sequential);CHKERRQ(ierr);
956   ierr = PetscFree(sub_schurs->auxglobal_parallel);CHKERRQ(ierr);
957   sub_schurs->n_subs = 0;
958   PetscFunctionReturn(0);
959 }
960 
961 #undef __FUNCT__
962 #define __FUNCT__ "PCBDDCAdjGetNextLayer_Private"
963 PETSC_STATIC_INLINE PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt* queue_tip,PetscInt n_prev,PetscBT touched,PetscInt* xadj,PetscInt* adjncy,PetscInt* n_added)
964 {
965   PetscInt       i,j,n;
966   PetscErrorCode ierr;
967 
968   PetscFunctionBegin;
969   n = 0;
970   for (i=-n_prev;i<0;i++) {
971     PetscInt start_dof = queue_tip[i];
972     for (j=xadj[start_dof];j<xadj[start_dof+1];j++) {
973       PetscInt dof = adjncy[j];
974       if (!PetscBTLookup(touched,dof)) {
975         ierr = PetscBTSet(touched,dof);CHKERRQ(ierr);
976         queue_tip[n] = dof;
977         n++;
978       }
979     }
980   }
981   *n_added = n;
982   PetscFunctionReturn(0);
983 }
984