xref: /petsc/src/ksp/pc/impls/bddc/bddcschurs.c (revision a64f4aa4610dd4ec03fac46179ae53980bc33f9b)
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 static PetscErrorCode PCBDDCMumpsInteriorDestroy(PC);
8 static PetscErrorCode PCBDDCMumpsInteriorSolve(PC,Vec,Vec);
9 
10 #undef __FUNCT__
11 #define __FUNCT__ "PCBDDCMumpsInteriorSolve"
12 static PetscErrorCode PCBDDCMumpsInteriorSolve(PC pc, Vec rhs, Vec sol)
13 {
14   PCBDDCMumpsInterior ctx;
15   PetscScalar         *array,*array_mumps;
16   PetscInt            ival;
17   PetscErrorCode      ierr;
18 
19   PetscFunctionBegin;
20   ierr = PCShellGetContext(pc,(void **)&ctx);CHKERRQ(ierr);
21   ierr = MatMumpsGetIcntl(ctx->F,26,&ival);CHKERRQ(ierr);
22   ierr = MatMumpsSetIcntl(ctx->F,26,0);CHKERRQ(ierr);
23   /* copy rhs into factored matrix workspace (can it be avoided?, MatSolve_MUMPS has another copy b->x internally) */
24   ierr = VecGetArrayRead(rhs,(const PetscScalar**)&array);CHKERRQ(ierr);
25   ierr = VecGetArray(ctx->rhs,&array_mumps);CHKERRQ(ierr);
26   ierr = PetscMemcpy(array_mumps,array,ctx->n*sizeof(PetscScalar));CHKERRQ(ierr);
27   ierr = VecRestoreArray(ctx->rhs,&array_mumps);CHKERRQ(ierr);
28   ierr = VecRestoreArrayRead(rhs,(const PetscScalar**)&array);CHKERRQ(ierr);
29 
30   ierr = MatSolve(ctx->F,ctx->rhs,ctx->sol);CHKERRQ(ierr);
31 
32   /* get back data to caller worskpace */
33   ierr = VecGetArrayRead(ctx->sol,(const PetscScalar**)&array_mumps);CHKERRQ(ierr);
34   ierr = VecGetArray(sol,&array);CHKERRQ(ierr);
35   ierr = PetscMemcpy(array,array_mumps,ctx->n*sizeof(PetscScalar));CHKERRQ(ierr);
36   ierr = VecRestoreArray(sol,&array);CHKERRQ(ierr);
37   ierr = VecRestoreArrayRead(ctx->sol,(const PetscScalar**)&array_mumps);CHKERRQ(ierr);
38   ierr = MatMumpsSetIcntl(ctx->F,26,ival);CHKERRQ(ierr);
39   PetscFunctionReturn(0);
40 }
41 
42 #undef __FUNCT__
43 #define __FUNCT__ "PCBDDCMumpsInteriorDestroy"
44 static PetscErrorCode PCBDDCMumpsInteriorDestroy(PC pc)
45 {
46   PCBDDCMumpsInterior ctx;
47   PetscErrorCode ierr;
48 
49   PetscFunctionBegin;
50   ierr = PCShellGetContext(pc,(void **)&ctx);CHKERRQ(ierr);
51   ierr = MatDestroy(&ctx->F);CHKERRQ(ierr);
52   ierr = VecDestroy(&ctx->sol);CHKERRQ(ierr);
53   ierr = VecDestroy(&ctx->rhs);CHKERRQ(ierr);
54   ierr = PetscFree(ctx);CHKERRQ(ierr);
55   PetscFunctionReturn(0);
56 }
57 
58 #undef __FUNCT__
59 #define __FUNCT__ "PCBDDCComputeExplicitSchur"
60 static PetscErrorCode PCBDDCComputeExplicitSchur(Mat M, PetscBool issym, MatReuse reuse, Mat *S)
61 {
62   Mat            B, C, D, Bd, Cd, AinvBd;
63   KSP            ksp;
64   PC             pc;
65   PetscBool      isLU, isILU, isCHOL, Bdense, Cdense;
66   PetscReal      fill = 2.0;
67   PetscInt       n_I;
68   PetscMPIInt    size;
69   PetscErrorCode ierr;
70 
71   PetscFunctionBegin;
72   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)M),&size);CHKERRQ(ierr);
73   if (size != 1) {
74     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not for parallel matrices");
75   }
76   if (reuse == MAT_REUSE_MATRIX) {
77     PetscBool Sdense;
78 
79     ierr = PetscObjectTypeCompare((PetscObject)*S, MATSEQDENSE, &Sdense);CHKERRQ(ierr);
80     if (!Sdense) {
81       SETERRQ(PetscObjectComm((PetscObject)M),PETSC_ERR_SUP,"S should dense");
82     }
83   }
84   ierr = MatSchurComplementGetSubMatrices(M, NULL, NULL, &B, &C, &D);CHKERRQ(ierr);
85   ierr = MatSchurComplementGetKSP(M, &ksp);CHKERRQ(ierr);
86   ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr);
87   ierr = PetscObjectTypeCompare((PetscObject) pc, PCLU, &isLU);CHKERRQ(ierr);
88   ierr = PetscObjectTypeCompare((PetscObject) pc, PCILU, &isILU);CHKERRQ(ierr);
89   ierr = PetscObjectTypeCompare((PetscObject) pc, PCCHOLESKY, &isCHOL);CHKERRQ(ierr);
90   ierr = PetscObjectTypeCompare((PetscObject) B, MATSEQDENSE, &Bdense);CHKERRQ(ierr);
91   ierr = PetscObjectTypeCompare((PetscObject) C, MATSEQDENSE, &Cdense);CHKERRQ(ierr);
92   ierr = MatGetSize(B,&n_I,NULL);CHKERRQ(ierr);
93   if (n_I) {
94     if (!Bdense) {
95       ierr = MatConvert(B, MATSEQDENSE, MAT_INITIAL_MATRIX, &Bd);CHKERRQ(ierr);
96     } else {
97       Bd = B;
98     }
99 
100     if (isLU || isILU || isCHOL) {
101       Mat fact;
102       ierr = KSPSetUp(ksp);CHKERRQ(ierr);
103       ierr = PCFactorGetMatrix(pc, &fact);CHKERRQ(ierr);
104       ierr = MatDuplicate(Bd, MAT_DO_NOT_COPY_VALUES, &AinvBd);CHKERRQ(ierr);
105       ierr = MatMatSolve(fact, Bd, AinvBd);CHKERRQ(ierr);
106     } else {
107       PetscBool ex = PETSC_TRUE;
108 
109       if (ex) {
110         Mat Ainvd;
111 
112         ierr = PCComputeExplicitOperator(pc, &Ainvd);CHKERRQ(ierr);
113         ierr = MatMatMult(Ainvd, Bd, MAT_INITIAL_MATRIX, fill, &AinvBd);CHKERRQ(ierr);
114         ierr = MatDestroy(&Ainvd);CHKERRQ(ierr);
115       } else {
116         Vec         sol,rhs;
117         PetscScalar *arrayrhs,*arraysol;
118         PetscInt    i,nrhs,n;
119 
120         ierr = MatDuplicate(Bd, MAT_DO_NOT_COPY_VALUES, &AinvBd);CHKERRQ(ierr);
121         ierr = MatGetSize(Bd,&n,&nrhs);CHKERRQ(ierr);
122         ierr = MatDenseGetArray(Bd,&arrayrhs);CHKERRQ(ierr);
123         ierr = MatDenseGetArray(AinvBd,&arraysol);CHKERRQ(ierr);
124         ierr = KSPGetSolution(ksp,&sol);CHKERRQ(ierr);
125         ierr = KSPGetRhs(ksp,&rhs);CHKERRQ(ierr);
126         for (i=0;i<nrhs;i++) {
127           ierr = VecPlaceArray(rhs,arrayrhs+i*n);CHKERRQ(ierr);
128           ierr = VecPlaceArray(sol,arraysol+i*n);CHKERRQ(ierr);
129           ierr = KSPSolve(ksp,rhs,sol);CHKERRQ(ierr);
130           ierr = VecResetArray(rhs);CHKERRQ(ierr);
131           ierr = VecResetArray(sol);CHKERRQ(ierr);
132         }
133         ierr = MatDenseRestoreArray(Bd,&arrayrhs);CHKERRQ(ierr);
134         ierr = MatDenseRestoreArray(AinvBd,&arrayrhs);CHKERRQ(ierr);
135       }
136     }
137     if (!Bdense & !issym) {
138       ierr = MatDestroy(&Bd);CHKERRQ(ierr);
139     }
140 
141     if (!issym) {
142       if (!Cdense) {
143         ierr = MatConvert(C, MATSEQDENSE, MAT_INITIAL_MATRIX, &Cd);CHKERRQ(ierr);
144       } else {
145         Cd = C;
146       }
147       ierr = MatMatMult(Cd, AinvBd, reuse, fill, S);CHKERRQ(ierr);
148       if (!Cdense) {
149         ierr = MatDestroy(&Cd);CHKERRQ(ierr);
150       }
151     } else {
152       ierr = MatTransposeMatMult(Bd, AinvBd, reuse, fill, S);CHKERRQ(ierr);
153       if (!Bdense) {
154         ierr = MatDestroy(&Bd);CHKERRQ(ierr);
155       }
156     }
157     ierr = MatDestroy(&AinvBd);CHKERRQ(ierr);
158   }
159 
160   if (D) {
161     Mat       Dd;
162     PetscBool Ddense;
163 
164     ierr = PetscObjectTypeCompare((PetscObject)D,MATSEQDENSE,&Ddense);CHKERRQ(ierr);
165     if (!Ddense) {
166       ierr = MatConvert(D, MATSEQDENSE, MAT_INITIAL_MATRIX, &Dd);CHKERRQ(ierr);
167     } else {
168       Dd = D;
169     }
170     if (n_I) {
171       ierr = MatAYPX(*S,-1.0,Dd,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
172     } else {
173       if (reuse == MAT_INITIAL_MATRIX) {
174         ierr = MatDuplicate(Dd,MAT_COPY_VALUES,S);CHKERRQ(ierr);
175       } else {
176         ierr = MatCopy(Dd,*S,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
177       }
178     }
179     if (!Ddense) {
180       ierr = MatDestroy(&Dd);CHKERRQ(ierr);
181     }
182   } else {
183     ierr = MatScale(*S,-1.0);CHKERRQ(ierr);
184   }
185   PetscFunctionReturn(0);
186 }
187 
188 #undef __FUNCT__
189 #define __FUNCT__ "PCBDDCSubSchursSetUp"
190 PetscErrorCode PCBDDCSubSchursSetUp(PCBDDCSubSchurs sub_schurs, Mat Ain, Mat Sin, PetscInt xadj[], PetscInt adjncy[], PetscInt nlayers, PetscBool faster_deluxe, PetscBool compute_Stilda, PetscBool use_edges, PetscBool use_faces)
191 {
192   Mat                    A_II,A_IB,A_BI,A_BB,AE_II;
193   Mat                    S_all,S_all_inv;
194   Mat                    global_schur_subsets,work_mat;
195   Mat                    S_Ej_tilda_all,S_Ej_inv_all;
196   ISLocalToGlobalMapping l2gmap_subsets;
197   IS                     is_I,temp_is;
198   PetscInt               *nnz,*all_local_idx_G,*all_local_idx_N;
199   PetscInt               *auxnum1,*auxnum2,*all_local_idx_G_rep;
200   PetscInt               i,subset_size,max_subset_size;
201   PetscInt               extra,local_size,global_size;
202   PetscBLASInt           B_N,B_ierr,B_lwork,*pivots;
203   PetscScalar            *Bwork;
204   PetscSubcomm           subcomm;
205   PetscMPIInt            color,rank;
206   MPI_Comm               comm_n;
207   PetscErrorCode         ierr;
208 
209   PetscFunctionBegin;
210   /* update info in sub_schurs */
211   ierr = MatDestroy(&sub_schurs->A);CHKERRQ(ierr);
212   ierr = MatDestroy(&sub_schurs->S);CHKERRQ(ierr);
213   if (Ain) {
214     PetscBool isseqaij;
215 
216     ierr = PetscObjectTypeCompare((PetscObject)Ain,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
217     if (isseqaij) {
218       ierr = PetscObjectReference((PetscObject)Ain);CHKERRQ(ierr);
219       sub_schurs->A = Ain;
220     } else { /* SeqBAIJ matrices does not support symmetry checking, SeqSBAIJ does not support MatPermute */
221       ierr = MatConvert(Ain,MATSEQAIJ,MAT_INITIAL_MATRIX,&sub_schurs->A);CHKERRQ(ierr);
222     }
223   }
224   ierr = PetscObjectReference((PetscObject)Sin);CHKERRQ(ierr);
225   sub_schurs->S = Sin;
226   if (sub_schurs->use_mumps) {
227     sub_schurs->use_mumps = (PetscBool)(!!sub_schurs->A);
228   }
229 
230   /* preliminary checks */
231   if (!sub_schurs->use_mumps && compute_Stilda) {
232     SETERRQ(PetscObjectComm((PetscObject)sub_schurs->l2gmap),PETSC_ERR_SUP,"Adaptive selection of constraints requires MUMPS");
233   }
234   /* determine if we are dealing with hermitian positive definite problems */
235   sub_schurs->is_hermitian = PETSC_FALSE;
236   sub_schurs->is_posdef = PETSC_FALSE;
237   if (sub_schurs->A) {
238     PetscInt lsize;
239 
240     ierr = MatGetSize(sub_schurs->A,&lsize,NULL);CHKERRQ(ierr);
241     if (lsize) {
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     } else {
257       sub_schurs->is_hermitian = PETSC_TRUE;
258       sub_schurs->is_posdef = PETSC_TRUE;
259     }
260     if (compute_Stilda && (!sub_schurs->is_hermitian || !sub_schurs->is_posdef)) {
261       SETERRQ(PetscObjectComm((PetscObject)sub_schurs->l2gmap),PETSC_ERR_SUP,"General matrix pencils are not currently supported");
262     }
263   }
264   /* restrict work on active processes */
265   color = 0;
266   if (!sub_schurs->n_subs) color = 1; /* this can happen if we are in a multilevel case or if the subdomain is disconnected */
267   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&rank);CHKERRQ(ierr);
268   ierr = PetscSubcommCreate(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&subcomm);CHKERRQ(ierr);
269   ierr = PetscSubcommSetNumber(subcomm,2);CHKERRQ(ierr);
270   ierr = PetscSubcommSetTypeGeneral(subcomm,color,rank);CHKERRQ(ierr);
271   ierr = PetscCommDuplicate(PetscSubcommChild(subcomm),&comm_n,NULL);CHKERRQ(ierr);
272   ierr = PetscSubcommDestroy(&subcomm);CHKERRQ(ierr);
273   if (!sub_schurs->n_subs) {
274     ierr = PetscCommDestroy(&comm_n);CHKERRQ(ierr);
275     PetscFunctionReturn(0);
276   }
277 
278   /* get Schur complement matrices */
279   if (!sub_schurs->use_mumps) {
280     Mat       tA_IB,tA_BI,tA_BB;
281     PetscBool isseqaij;
282     ierr = MatSchurComplementGetSubMatrices(sub_schurs->S,&A_II,NULL,&tA_IB,&tA_BI,&tA_BB);CHKERRQ(ierr);
283     ierr = PetscObjectTypeCompare((PetscObject)tA_BB,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
284     if (!isseqaij) {
285       ierr = MatConvert(tA_BB,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_BB);CHKERRQ(ierr);
286       ierr = MatConvert(tA_IB,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_IB);CHKERRQ(ierr);
287       ierr = MatConvert(tA_BI,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_BI);CHKERRQ(ierr);
288     } else {
289       ierr = PetscObjectReference((PetscObject)tA_BB);CHKERRQ(ierr);
290       A_BB = tA_BB;
291       ierr = PetscObjectReference((PetscObject)tA_IB);CHKERRQ(ierr);
292       A_IB = tA_IB;
293       ierr = PetscObjectReference((PetscObject)tA_BI);CHKERRQ(ierr);
294       A_BI = tA_BI;
295     }
296   } else {
297     A_II = NULL;
298     A_IB = NULL;
299     A_BI = NULL;
300     A_BB = NULL;
301   }
302   S_all = NULL;
303   S_all_inv = NULL;
304   S_Ej_tilda_all = NULL;
305   S_Ej_inv_all = NULL;
306 
307   /* determine interior problems */
308   ierr = ISDestroy(&sub_schurs->is_I_layer);CHKERRQ(ierr);
309   ierr = ISGetLocalSize(sub_schurs->is_I,&i);CHKERRQ(ierr);
310   if (nlayers >= 0 && i) { /* Interior problems can be different from the original one */
311     PetscBT                touched;
312     const PetscInt*        idx_B;
313     PetscInt               n_I,n_B,n_local_dofs,n_prev_added,j,layer,*local_numbering;
314 
315     if (xadj == NULL || adjncy == NULL) {
316       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cannot request layering without adjacency");
317     }
318     /* get sizes */
319     ierr = ISGetLocalSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr);
320     ierr = ISGetLocalSize(sub_schurs->is_B,&n_B);CHKERRQ(ierr);
321 
322     ierr = PetscMalloc1(n_I+n_B,&local_numbering);CHKERRQ(ierr);
323     ierr = PetscBTCreate(n_I+n_B,&touched);CHKERRQ(ierr);
324     ierr = PetscBTMemzero(n_I+n_B,touched);CHKERRQ(ierr);
325 
326     /* all boundary dofs must be skipped when adding layers */
327     ierr = ISGetIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr);
328     for (j=0;j<n_B;j++) {
329       ierr = PetscBTSet(touched,idx_B[j]);CHKERRQ(ierr);
330     }
331     ierr = PetscMemcpy(local_numbering,idx_B,n_B*sizeof(PetscInt));CHKERRQ(ierr);
332     ierr = ISRestoreIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr);
333 
334     /* add prescribed number of layers of dofs */
335     n_local_dofs = n_B;
336     n_prev_added = n_B;
337     for (layer=0;layer<nlayers;layer++) {
338       PetscInt n_added;
339       if (n_local_dofs == n_I+n_B) break;
340       if (n_local_dofs > n_I+n_B) {
341         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);
342       }
343       ierr = PCBDDCAdjGetNextLayer_Private(local_numbering+n_local_dofs,n_prev_added,touched,xadj,adjncy,&n_added);CHKERRQ(ierr);
344       n_prev_added = n_added;
345       n_local_dofs += n_added;
346       if (!n_added) break;
347     }
348     ierr = PetscBTDestroy(&touched);CHKERRQ(ierr);
349 
350     /* IS for I layer dofs in original numbering */
351     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);
352     ierr = PetscFree(local_numbering);CHKERRQ(ierr);
353     ierr = ISSort(sub_schurs->is_I_layer);CHKERRQ(ierr);
354     /* IS for I layer dofs in I numbering */
355     if (!sub_schurs->use_mumps) {
356       ISLocalToGlobalMapping ItoNmap;
357       ierr = ISLocalToGlobalMappingCreateIS(sub_schurs->is_I,&ItoNmap);CHKERRQ(ierr);
358       ierr = ISGlobalToLocalMappingApplyIS(ItoNmap,IS_GTOLM_DROP,sub_schurs->is_I_layer,&is_I);CHKERRQ(ierr);
359       ierr = ISLocalToGlobalMappingDestroy(&ItoNmap);CHKERRQ(ierr);
360 
361       /* II block */
362       ierr = MatGetSubMatrix(A_II,is_I,is_I,MAT_INITIAL_MATRIX,&AE_II);CHKERRQ(ierr);
363     }
364   } else {
365     PetscInt n_I;
366 
367     /* IS for I dofs in original numbering */
368     ierr = PetscObjectReference((PetscObject)sub_schurs->is_I);CHKERRQ(ierr);
369     sub_schurs->is_I_layer = sub_schurs->is_I;
370 
371     /* IS for I dofs in I numbering (strided 1) */
372     if (!sub_schurs->use_mumps) {
373       ierr = ISGetSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr);
374       ierr = ISCreateStride(PetscObjectComm((PetscObject)sub_schurs->is_I),n_I,0,1,&is_I);CHKERRQ(ierr);
375 
376       /* II block is the same */
377       ierr = PetscObjectReference((PetscObject)A_II);CHKERRQ(ierr);
378       AE_II = A_II;
379     }
380   }
381 
382   /* Get info on subset sizes and sum of all subsets sizes */
383   max_subset_size = 0;
384   local_size = 0;
385   for (i=0;i<sub_schurs->n_subs;i++) {
386     ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
387     max_subset_size = PetscMax(subset_size,max_subset_size);
388     local_size += subset_size;
389   }
390 
391   /* Work arrays for local indices */
392   extra = 0;
393   if (sub_schurs->use_mumps) {
394     ierr = ISGetLocalSize(sub_schurs->is_I_layer,&extra);CHKERRQ(ierr);
395   }
396   ierr = PetscMalloc1(local_size+extra,&all_local_idx_N);CHKERRQ(ierr);
397   if (extra) {
398     const PetscInt *idxs;
399     ierr = ISGetIndices(sub_schurs->is_I_layer,&idxs);CHKERRQ(ierr);
400     ierr = PetscMemcpy(all_local_idx_N,idxs,extra*sizeof(PetscInt));CHKERRQ(ierr);
401     ierr = ISRestoreIndices(sub_schurs->is_I_layer,&idxs);CHKERRQ(ierr);
402   }
403   ierr = PetscMalloc1(local_size,&nnz);CHKERRQ(ierr);
404   ierr = PetscMalloc2(sub_schurs->n_subs,&auxnum1,sub_schurs->n_subs,&auxnum2);CHKERRQ(ierr);
405 
406   /* Get local indices in local numbering */
407   local_size = 0;
408   for (i=0;i<sub_schurs->n_subs;i++) {
409     PetscInt j;
410     const    PetscInt *idxs;
411 
412     ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
413     ierr = ISGetIndices(sub_schurs->is_subs[i],&idxs);CHKERRQ(ierr);
414     /* start (smallest in global ordering) and multiplicity */
415     auxnum1[i] = idxs[0];
416     auxnum2[i] = subset_size;
417     /* subset indices in local numbering */
418     ierr = PetscMemcpy(all_local_idx_N+local_size+extra,idxs,subset_size*sizeof(PetscInt));CHKERRQ(ierr);
419     ierr = ISRestoreIndices(sub_schurs->is_subs[i],&idxs);CHKERRQ(ierr);
420     for (j=0;j<subset_size;j++) nnz[local_size+j] = subset_size;
421     local_size += subset_size;
422   }
423 
424   /* allocate extra workspace needed only for GETRI */
425   Bwork = NULL;
426   pivots = NULL;
427   if (sub_schurs->n_subs && !sub_schurs->is_hermitian) {
428     PetscScalar lwork;
429 
430     B_lwork = -1;
431     ierr = PetscBLASIntCast(local_size,&B_N);CHKERRQ(ierr);
432     ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
433     PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,Bwork,&B_N,pivots,&lwork,&B_lwork,&B_ierr));
434     ierr = PetscFPTrapPop();CHKERRQ(ierr);
435     if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to GETRI Lapack routine %d",(int)B_ierr);
436     ierr = PetscBLASIntCast((PetscInt)PetscRealPart(lwork),&B_lwork);CHKERRQ(ierr);
437     ierr = PetscMalloc2(B_lwork,&Bwork,B_N,&pivots);CHKERRQ(ierr);
438   }
439 
440   /* prepare parallel matrices for summing up properly schurs on subsets */
441   ierr = PCBDDCSubsetNumbering(comm_n,sub_schurs->l2gmap,sub_schurs->n_subs,auxnum1,auxnum2,&global_size,&all_local_idx_G_rep);CHKERRQ(ierr);
442   ierr = PetscMalloc1(local_size,&all_local_idx_G);CHKERRQ(ierr);
443   local_size = 0;
444   for (i=0;i<sub_schurs->n_subs;i++) {
445     PetscInt j;
446     for (j=0;j<auxnum2[i];j++) all_local_idx_G[local_size++] = all_local_idx_G_rep[i] + j;
447   }
448   ierr = PetscFree(all_local_idx_G_rep);CHKERRQ(ierr);
449   ierr = PetscFree2(auxnum1,auxnum2);CHKERRQ(ierr);
450   ierr = ISLocalToGlobalMappingCreate(comm_n,1,local_size,all_local_idx_G,PETSC_COPY_VALUES,&l2gmap_subsets);CHKERRQ(ierr);
451   ierr = MatCreateIS(comm_n,1,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size,l2gmap_subsets,&work_mat);CHKERRQ(ierr);
452   ierr = ISLocalToGlobalMappingDestroy(&l2gmap_subsets);CHKERRQ(ierr);
453   ierr = MatCreate(PetscObjectComm((PetscObject)work_mat),&global_schur_subsets);CHKERRQ(ierr);
454   ierr = MatSetSizes(global_schur_subsets,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size);CHKERRQ(ierr);
455   ierr = MatSetType(global_schur_subsets,MATMPIAIJ);CHKERRQ(ierr);
456 
457   /* subset indices in local boundary numbering */
458   if (!sub_schurs->is_Ej_all) {
459     PetscInt *all_local_idx_B;
460 
461     ierr = PetscMalloc1(local_size,&all_local_idx_B);CHKERRQ(ierr);
462     ierr = ISGlobalToLocalMappingApply(sub_schurs->BtoNmap,IS_GTOLM_DROP,local_size,all_local_idx_N+extra,&subset_size,all_local_idx_B);CHKERRQ(ierr);
463     if (subset_size != local_size) {
464       SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in sub_schurs serial (BtoNmap)! %d != %d\n",subset_size,local_size);
465     }
466     ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size,all_local_idx_B,PETSC_OWN_POINTER,&sub_schurs->is_Ej_all);CHKERRQ(ierr);
467   }
468 
469   /* Local matrix of all local Schur on subsets (transposed) */
470   if (!sub_schurs->S_Ej_all) {
471     ierr = MatCreate(PETSC_COMM_SELF,&sub_schurs->S_Ej_all);CHKERRQ(ierr);
472     ierr = MatSetSizes(sub_schurs->S_Ej_all,PETSC_DECIDE,PETSC_DECIDE,local_size,local_size);CHKERRQ(ierr);
473     ierr = MatSetType(sub_schurs->S_Ej_all,MATAIJ);CHKERRQ(ierr);
474     ierr = MatSeqAIJSetPreallocation(sub_schurs->S_Ej_all,0,nnz);CHKERRQ(ierr);
475   } else {
476     ierr = MatZeroEntries(sub_schurs->S_Ej_all);CHKERRQ(ierr);
477   }
478 
479   /* Compute Schur complements explicitly */
480   ierr = PetscBTMemzero(sub_schurs->n_subs,sub_schurs->computed_Stilda_subs);CHKERRQ(ierr);
481   if (!sub_schurs->use_mumps) {
482     Mat         S_Ej_expl;
483     PetscScalar *work;
484     PetscInt    j,*dummy_idx;
485     PetscBool   Sdense;
486 
487     ierr = PetscMalloc2(max_subset_size,&dummy_idx,max_subset_size*max_subset_size,&work);CHKERRQ(ierr);
488     local_size = 0;
489     for (i=0;i<sub_schurs->n_subs;i++) {
490       IS  is_subset_B;
491       Mat AE_EE,AE_IE,AE_EI,S_Ej;
492 
493       /* subsets in original and boundary numbering */
494       ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,sub_schurs->is_subs[i],&is_subset_B);CHKERRQ(ierr);
495       /* EE block */
496       ierr = MatGetSubMatrix(A_BB,is_subset_B,is_subset_B,MAT_INITIAL_MATRIX,&AE_EE);CHKERRQ(ierr);
497       /* IE block */
498       ierr = MatGetSubMatrix(A_IB,is_I,is_subset_B,MAT_INITIAL_MATRIX,&AE_IE);CHKERRQ(ierr);
499       /* EI block */
500       if (sub_schurs->is_hermitian) {
501         ierr = MatCreateTranspose(AE_IE,&AE_EI);CHKERRQ(ierr);
502       } else {
503         ierr = MatGetSubMatrix(A_BI,is_subset_B,is_I,MAT_INITIAL_MATRIX,&AE_EI);CHKERRQ(ierr);
504       }
505       ierr = ISDestroy(&is_subset_B);CHKERRQ(ierr);
506       ierr = MatCreateSchurComplement(AE_II,AE_II,AE_IE,AE_EI,AE_EE,&S_Ej);CHKERRQ(ierr);
507       ierr = MatDestroy(&AE_EE);CHKERRQ(ierr);
508       ierr = MatDestroy(&AE_IE);CHKERRQ(ierr);
509       ierr = MatDestroy(&AE_EI);CHKERRQ(ierr);
510       if (AE_II == A_II) { /* we can reuse the same ksp */
511         KSP ksp;
512         ierr = MatSchurComplementGetKSP(sub_schurs->S,&ksp);CHKERRQ(ierr);
513         ierr = MatSchurComplementSetKSP(S_Ej,ksp);CHKERRQ(ierr);
514       } else { /* build new ksp object which inherits ksp and pc types from the original one */
515         KSP       origksp,schurksp;
516         PC        origpc,schurpc;
517         KSPType   ksp_type;
518         PetscInt  n_internal;
519         PetscBool ispcnone;
520 
521         ierr = MatSchurComplementGetKSP(sub_schurs->S,&origksp);CHKERRQ(ierr);
522         ierr = MatSchurComplementGetKSP(S_Ej,&schurksp);CHKERRQ(ierr);
523         ierr = KSPGetType(origksp,&ksp_type);CHKERRQ(ierr);
524         ierr = KSPSetType(schurksp,ksp_type);CHKERRQ(ierr);
525         ierr = KSPGetPC(schurksp,&schurpc);CHKERRQ(ierr);
526         ierr = KSPGetPC(origksp,&origpc);CHKERRQ(ierr);
527         ierr = PetscObjectTypeCompare((PetscObject)origpc,PCNONE,&ispcnone);CHKERRQ(ierr);
528         if (!ispcnone) {
529           PCType pc_type;
530           ierr = PCGetType(origpc,&pc_type);CHKERRQ(ierr);
531           ierr = PCSetType(schurpc,pc_type);CHKERRQ(ierr);
532         } else {
533           ierr = PCSetType(schurpc,PCLU);CHKERRQ(ierr);
534         }
535         ierr = ISGetSize(is_I,&n_internal);CHKERRQ(ierr);
536         if (n_internal) { /* UMFPACK gives error with 0 sized problems */
537           MatSolverPackage solver=NULL;
538           ierr = PCFactorGetMatSolverPackage(origpc,(const MatSolverPackage*)&solver);CHKERRQ(ierr);
539           if (solver) {
540             ierr = PCFactorSetMatSolverPackage(schurpc,solver);CHKERRQ(ierr);
541           }
542         }
543         ierr = KSPSetUp(schurksp);CHKERRQ(ierr);
544       }
545       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
546       ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&S_Ej_expl);CHKERRQ(ierr);
547       ierr = PCBDDCComputeExplicitSchur(S_Ej,sub_schurs->is_hermitian,MAT_REUSE_MATRIX,&S_Ej_expl);CHKERRQ(ierr);
548       ierr = PetscObjectTypeCompare((PetscObject)S_Ej_expl,MATSEQDENSE,&Sdense);CHKERRQ(ierr);
549       if (Sdense) {
550         for (j=0;j<subset_size;j++) {
551           dummy_idx[j]=local_size+j;
552         }
553         ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
554       } else {
555         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not yet implemented for sparse matrices");
556       }
557       ierr = MatDestroy(&S_Ej);CHKERRQ(ierr);
558       ierr = MatDestroy(&S_Ej_expl);CHKERRQ(ierr);
559       local_size += subset_size;
560     }
561     ierr = PetscFree2(dummy_idx,work);CHKERRQ(ierr);
562     /* free */
563     ierr = ISDestroy(&is_I);CHKERRQ(ierr);
564     ierr = MatDestroy(&AE_II);CHKERRQ(ierr);
565     ierr = PetscFree(all_local_idx_N);CHKERRQ(ierr);
566   } else {
567     Mat         A,F;
568     IS          is_A_all;
569     PetscScalar *work;
570     PetscInt    *idxs_schur,n_I,n_I_all,*dummy_idx;
571 
572     /* get working mat */
573     ierr = ISGetLocalSize(sub_schurs->is_I_layer,&n_I);CHKERRQ(ierr);
574     ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size+n_I,all_local_idx_N,PETSC_COPY_VALUES,&is_A_all);CHKERRQ(ierr);
575     ierr = MatGetSubMatrixUnsorted(sub_schurs->A,is_A_all,is_A_all,&A);CHKERRQ(ierr);
576     ierr = ISDestroy(&is_A_all);CHKERRQ(ierr);
577 
578     if (n_I) {
579       if (sub_schurs->is_hermitian && sub_schurs->is_posdef) {
580         ierr = MatGetFactor(A,MATSOLVERMUMPS,MAT_FACTOR_CHOLESKY,&F);CHKERRQ(ierr);
581       } else {
582         ierr = MatGetFactor(A,MATSOLVERMUMPS,MAT_FACTOR_LU,&F);CHKERRQ(ierr);
583       }
584 
585       /* subsets ordered last */
586       ierr = PetscMalloc1(local_size,&idxs_schur);CHKERRQ(ierr);
587       for (i=0;i<local_size;i++) {
588         idxs_schur[i] = n_I+i+1;
589       }
590 #if defined(PETSC_HAVE_MUMPS)
591       ierr = MatMumpsSetSchurIndices(F,local_size,idxs_schur);CHKERRQ(ierr);
592 #endif
593       ierr = PetscFree(idxs_schur);CHKERRQ(ierr);
594 
595       /* factorization step */
596       if (sub_schurs->is_hermitian && sub_schurs->is_posdef) {
597         ierr = MatCholeskyFactorSymbolic(F,A,NULL,NULL);CHKERRQ(ierr);
598         ierr = MatCholeskyFactorNumeric(F,A,NULL);CHKERRQ(ierr);
599       } else {
600         ierr = MatLUFactorSymbolic(F,A,NULL,NULL,NULL);CHKERRQ(ierr);
601         ierr = MatLUFactorNumeric(F,A,NULL);CHKERRQ(ierr);
602       }
603 
604       /* get explicit Schur Complement computed during numeric factorization */
605 #if defined(PETSC_HAVE_MUMPS)
606       ierr = MatMumpsGetSchurComplement(F,&S_all);CHKERRQ(ierr);
607 #endif
608 
609       /* we can reuse the interior solver if we are not using the economic version */
610       ierr = ISGetLocalSize(sub_schurs->is_I,&n_I_all);CHKERRQ(ierr);
611       if (n_I == n_I_all) {
612         PCBDDCMumpsInterior msolv_ctx;
613 
614         ierr = PCDestroy(&sub_schurs->interior_solver);CHKERRQ(ierr);
615         ierr = PetscNew(&msolv_ctx);CHKERRQ(ierr);
616         msolv_ctx->n = n_I;
617         ierr = PetscObjectReference((PetscObject)F);CHKERRQ(ierr);
618         msolv_ctx->F = F;
619         ierr = MatCreateVecs(F,&msolv_ctx->sol,&msolv_ctx->rhs);CHKERRQ(ierr);
620         ierr = PCCreate(PETSC_COMM_SELF,&sub_schurs->interior_solver);CHKERRQ(ierr);
621         ierr = MatSchurComplementGetSubMatrices(sub_schurs->S,&A_II,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
622         ierr = PCSetOperators(sub_schurs->interior_solver,A_II,A_II);CHKERRQ(ierr);
623         ierr = PCSetType(sub_schurs->interior_solver,PCSHELL);CHKERRQ(ierr);
624         ierr = PCShellSetContext(sub_schurs->interior_solver,msolv_ctx);CHKERRQ(ierr);
625         ierr = PCShellSetApply(sub_schurs->interior_solver,PCBDDCMumpsInteriorSolve);CHKERRQ(ierr);
626         ierr = PCShellSetDestroy(sub_schurs->interior_solver,PCBDDCMumpsInteriorDestroy);CHKERRQ(ierr);
627       }
628       ierr = MatDestroy(&F);CHKERRQ(ierr);
629     } else {
630       ierr = MatConvert(A,MATSEQDENSE,MAT_INITIAL_MATRIX,&S_all);CHKERRQ(ierr);
631     }
632     ierr = MatDestroy(&A);CHKERRQ(ierr);
633     ierr = PetscFree(all_local_idx_N);CHKERRQ(ierr);
634 
635     if (compute_Stilda) {
636       ierr = MatCreate(PETSC_COMM_SELF,&S_Ej_tilda_all);CHKERRQ(ierr);
637       ierr = MatSetSizes(S_Ej_tilda_all,PETSC_DECIDE,PETSC_DECIDE,local_size,local_size);CHKERRQ(ierr);
638       ierr = MatSetType(S_Ej_tilda_all,MATAIJ);CHKERRQ(ierr);
639       ierr = MatSeqAIJSetPreallocation(S_Ej_tilda_all,0,nnz);CHKERRQ(ierr);
640       ierr = MatCreate(PETSC_COMM_SELF,&S_Ej_inv_all);CHKERRQ(ierr);
641       ierr = MatSetSizes(S_Ej_inv_all,PETSC_DECIDE,PETSC_DECIDE,local_size,local_size);CHKERRQ(ierr);
642       ierr = MatSetType(S_Ej_inv_all,MATAIJ);CHKERRQ(ierr);
643       ierr = MatSeqAIJSetPreallocation(S_Ej_inv_all,0,nnz);CHKERRQ(ierr);
644 
645       /* compute St^-1 */
646       if (local_size) { /* multilevel guard */
647         PetscScalar *vals;
648 
649         ierr = PetscBLASIntCast(local_size,&B_N);CHKERRQ(ierr);
650         ierr = MatDuplicate(S_all,MAT_COPY_VALUES,&S_all_inv);CHKERRQ(ierr);
651         ierr = MatDenseGetArray(S_all_inv,&vals);CHKERRQ(ierr);
652         if (!sub_schurs->is_hermitian) {
653           PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,vals,&B_N,pivots,&B_ierr));
654           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
655           PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,vals,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
656           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
657         } else {
658           PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,vals,&B_N,&B_ierr));
659           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
660           PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,vals,&B_N,&B_ierr));
661           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
662         }
663         ierr = MatDenseRestoreArray(S_all_inv,&vals);CHKERRQ(ierr);
664       }
665     }
666 
667     /* Work arrays */
668     if (sub_schurs->n_subs == 1) {
669       ierr = PetscMalloc1(max_subset_size,&dummy_idx);CHKERRQ(ierr);
670     } else {
671       ierr = PetscMalloc2(max_subset_size,&dummy_idx,max_subset_size*max_subset_size,&work);CHKERRQ(ierr);
672     }
673 
674     local_size = 0;
675     for (i=0;i<sub_schurs->n_subs;i++) {
676       Mat S_Ej;
677       IS  is_E;
678       PetscInt j;
679 
680       /* get S_E */
681       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
682       if (sub_schurs->n_subs == 1) {
683         ierr = MatDenseGetArray(S_all,&work);CHKERRQ(ierr);
684         S_Ej = NULL;
685         is_E = NULL;
686       } else {
687         ierr = ISCreateStride(PETSC_COMM_SELF,subset_size,local_size,1,&is_E);CHKERRQ(ierr);
688         ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&S_Ej);CHKERRQ(ierr);
689         ierr = MatGetSubMatrix(S_all,is_E,is_E,MAT_REUSE_MATRIX,&S_Ej);CHKERRQ(ierr);
690       }
691       /* insert S_E values */
692       for (j=0;j<subset_size;j++) {
693         dummy_idx[j]=local_size+j;
694       }
695       ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
696 
697       /* if adaptivity is requested, invert S_E and insert St_E^-1 blocks */
698       if (compute_Stilda && ((PetscBTLookup(sub_schurs->is_edge,i) && use_edges) || (!PetscBTLookup(sub_schurs->is_edge,i) && use_faces))) {
699         /* get S_E^-1 */
700         ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr);
701         ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
702         if (!sub_schurs->is_hermitian) {
703           PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,work,&B_N,pivots,&B_ierr));
704           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
705           PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,work,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
706           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
707         } else {
708           PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,work,&B_N,&B_ierr));
709           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
710           PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,work,&B_N,&B_ierr));
711           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
712         }
713         ierr = PetscFPTrapPop();CHKERRQ(ierr);
714         ierr = MatSetValues(S_Ej_inv_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
715 
716         /* get St_E^-1 */
717         if (sub_schurs->n_subs == 1) {
718           ierr = MatDenseRestoreArray(S_all,&work);CHKERRQ(ierr);
719           ierr = MatDenseGetArray(S_all_inv,&work);CHKERRQ(ierr);
720         } else {
721           ierr = MatGetSubMatrix(S_all_inv,is_E,is_E,MAT_REUSE_MATRIX,&S_Ej);CHKERRQ(ierr);
722         }
723         ierr = MatSetValues(S_Ej_tilda_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
724         if (sub_schurs->n_subs == 1) {
725           ierr = MatDenseRestoreArray(S_all_inv,&work);CHKERRQ(ierr);
726         }
727         ierr = PetscBTSet(sub_schurs->computed_Stilda_subs,i);CHKERRQ(ierr);
728       } else if (sub_schurs->n_subs == 1) {
729         ierr = MatDenseRestoreArray(S_all,&work);CHKERRQ(ierr);
730       }
731       ierr = MatDestroy(&S_Ej);CHKERRQ(ierr);
732       ierr = ISDestroy(&is_E);CHKERRQ(ierr);
733       local_size += subset_size;
734     }
735     if (sub_schurs->n_subs == 1) {
736       ierr = PetscFree(dummy_idx);CHKERRQ(ierr);
737     } else {
738       ierr = PetscFree2(dummy_idx,work);CHKERRQ(ierr);
739     }
740   }
741   ierr = PetscFree(nnz);CHKERRQ(ierr);
742   ierr = MatDestroy(&S_all);CHKERRQ(ierr);
743   ierr = MatDestroy(&S_all_inv);CHKERRQ(ierr);
744   ierr = MatDestroy(&A_BB);CHKERRQ(ierr);
745   ierr = MatDestroy(&A_IB);CHKERRQ(ierr);
746   ierr = MatDestroy(&A_BI);CHKERRQ(ierr);
747   ierr = MatAssemblyBegin(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
748   ierr = MatAssemblyEnd(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
749   if (compute_Stilda) {
750     ierr = MatAssemblyBegin(S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
751     ierr = MatAssemblyEnd(S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
752     ierr = MatAssemblyBegin(S_Ej_inv_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
753     ierr = MatAssemblyEnd(S_Ej_inv_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
754   }
755 
756   /* Global matrix of all assembled Schur on subsets */
757   ierr = MatISSetLocalMat(work_mat,sub_schurs->S_Ej_all);CHKERRQ(ierr);
758   ierr = MatISSetMPIXAIJPreallocation_Private(work_mat,global_schur_subsets,PETSC_TRUE);CHKERRQ(ierr);
759   ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
760 
761   /* Get local part of (\sum_j S_Ej) */
762   ierr = ISCreateGeneral(comm_n,local_size,all_local_idx_G,PETSC_OWN_POINTER,&temp_is);CHKERRQ(ierr);
763   ierr = MatDestroy(&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
764   ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
765 
766   /* Compute explicitly (\sum_j S_Ej)^-1 (faster scaling during PCApply, needs extra work when doing setup) */
767   if (faster_deluxe) {
768     Mat         tmpmat;
769     PetscScalar *array;
770     PetscInt    cum;
771 
772     ierr = MatSeqAIJGetArray(sub_schurs->sum_S_Ej_all,&array);CHKERRQ(ierr);
773     cum = 0;
774     for (i=0;i<sub_schurs->n_subs;i++) {
775       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
776       ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr);
777       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
778       if (!sub_schurs->is_hermitian) {
779         PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,array+cum,&B_N,pivots,&B_ierr));
780         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
781         PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,array+cum,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
782         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
783       } else {
784         PetscInt j,k;
785 
786         PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,array+cum,&B_N,&B_ierr));
787         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
788         PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,array+cum,&B_N,&B_ierr));
789         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
790         for (j=0;j<B_N;j++) {
791           for (k=j+1;k<B_N;k++) {
792             array[k*B_N+j+cum] = array[j*B_N+k+cum];
793           }
794         }
795       }
796       ierr = PetscFPTrapPop();CHKERRQ(ierr);
797       cum += subset_size*subset_size;
798     }
799     ierr = MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_all,&array);CHKERRQ(ierr);
800     ierr = MatMatMult(sub_schurs->S_Ej_all,sub_schurs->sum_S_Ej_all,MAT_INITIAL_MATRIX,1.0,&tmpmat);CHKERRQ(ierr);
801     ierr = MatDestroy(&sub_schurs->S_Ej_all);CHKERRQ(ierr);
802     ierr = MatDestroy(&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
803     sub_schurs->S_Ej_all = tmpmat;
804   }
805 
806   /* Get local part of (\sum_j S^-1_Ej) (\sum_j St^-1_Ej) */
807   if (compute_Stilda) {
808     ierr = MatISSetLocalMat(work_mat,S_Ej_tilda_all);CHKERRQ(ierr);
809     ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
810     ierr = MatDestroy(&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
811     ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
812     ierr = MatISSetLocalMat(work_mat,S_Ej_inv_all);CHKERRQ(ierr);
813     ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
814     ierr = MatDestroy(&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
815     ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
816   }
817 
818   /* free workspace */
819   ierr = PetscFree2(Bwork,pivots);CHKERRQ(ierr);
820   ierr = MatDestroy(&global_schur_subsets);CHKERRQ(ierr);
821   ierr = MatDestroy(&S_Ej_tilda_all);CHKERRQ(ierr);
822   ierr = MatDestroy(&S_Ej_inv_all);CHKERRQ(ierr);
823   ierr = MatDestroy(&work_mat);CHKERRQ(ierr);
824   ierr = ISDestroy(&temp_is);CHKERRQ(ierr);
825   ierr = PetscCommDestroy(&comm_n);CHKERRQ(ierr);
826   PetscFunctionReturn(0);
827 }
828 
829 #undef __FUNCT__
830 #define __FUNCT__ "PCBDDCSubSchursInit"
831 PetscErrorCode PCBDDCSubSchursInit(PCBDDCSubSchurs sub_schurs, IS is_I, IS is_B, PCBDDCGraph graph, ISLocalToGlobalMapping BtoNmap)
832 {
833   IS              *faces,*edges,*all_cc,vertices;
834   PetscInt        i,n_faces,n_edges,n_all_cc;
835   PetscBool       is_sorted;
836   PetscErrorCode  ierr;
837 
838   PetscFunctionBegin;
839   ierr = ISSorted(is_I,&is_sorted);CHKERRQ(ierr);
840   if (!is_sorted) {
841     SETERRQ(PetscObjectComm((PetscObject)is_I),PETSC_ERR_PLIB,"IS for I dofs should be shorted");
842   }
843   ierr = ISSorted(is_B,&is_sorted);CHKERRQ(ierr);
844   if (!is_sorted) {
845     SETERRQ(PetscObjectComm((PetscObject)is_B),PETSC_ERR_PLIB,"IS for B dofs should be shorted");
846   }
847 
848   /* reset any previous data */
849   ierr = PCBDDCSubSchursReset(sub_schurs);CHKERRQ(ierr);
850 
851   /* get index sets for faces and edges (already sorted by global ordering) */
852   ierr = PCBDDCGraphGetCandidatesIS(graph,&n_faces,&faces,&n_edges,&edges,&vertices);CHKERRQ(ierr);
853   n_all_cc = n_faces+n_edges;
854   ierr = PetscBTCreate(n_all_cc,&sub_schurs->is_edge);CHKERRQ(ierr);
855   ierr = PetscBTCreate(n_all_cc,&sub_schurs->computed_Stilda_subs);CHKERRQ(ierr);
856   ierr = PetscMalloc1(n_all_cc,&all_cc);CHKERRQ(ierr);
857   for (i=0;i<n_faces;i++) {
858     all_cc[i] = faces[i];
859   }
860   for (i=0;i<n_edges;i++) {
861     all_cc[n_faces+i] = edges[i];
862     ierr = PetscBTSet(sub_schurs->is_edge,n_faces+i);CHKERRQ(ierr);
863   }
864   ierr = PetscFree(faces);CHKERRQ(ierr);
865   ierr = PetscFree(edges);CHKERRQ(ierr);
866 
867   /* Determine if MUMPS cane be used */
868   sub_schurs->use_mumps = PETSC_FALSE;
869 #if defined(PETSC_HAVE_MUMPS)
870   sub_schurs->use_mumps = PETSC_TRUE;
871 #endif
872 
873   ierr = PetscObjectReference((PetscObject)is_I);CHKERRQ(ierr);
874   sub_schurs->is_I = is_I;
875   ierr = PetscObjectReference((PetscObject)is_B);CHKERRQ(ierr);
876   sub_schurs->is_B = is_B;
877   ierr = PetscObjectReference((PetscObject)graph->l2gmap);CHKERRQ(ierr);
878   sub_schurs->l2gmap = graph->l2gmap;
879   ierr = PetscObjectReference((PetscObject)BtoNmap);CHKERRQ(ierr);
880   sub_schurs->BtoNmap = BtoNmap;
881   sub_schurs->n_subs = n_all_cc;
882   sub_schurs->is_subs = all_cc;
883   if (!sub_schurs->use_mumps) { /* sort by local ordering mumps is not present */
884     for (i=0;i<sub_schurs->n_subs;i++) {
885       ierr = ISSort(sub_schurs->is_subs[i]);CHKERRQ(ierr);
886     }
887   }
888   sub_schurs->is_Ej_com = vertices;
889 
890   /* allocate space for schur complements */
891   sub_schurs->S_Ej_all = NULL;
892   sub_schurs->sum_S_Ej_all = NULL;
893   sub_schurs->sum_S_Ej_inv_all = NULL;
894   sub_schurs->sum_S_Ej_tilda_all = NULL;
895   sub_schurs->is_Ej_all = NULL;
896   PetscFunctionReturn(0);
897 }
898 
899 #undef __FUNCT__
900 #define __FUNCT__ "PCBDDCSubSchursCreate"
901 PetscErrorCode PCBDDCSubSchursCreate(PCBDDCSubSchurs *sub_schurs)
902 {
903   PCBDDCSubSchurs schurs_ctx;
904   PetscErrorCode  ierr;
905 
906   PetscFunctionBegin;
907   ierr = PetscNew(&schurs_ctx);CHKERRQ(ierr);
908   schurs_ctx->n_subs = 0;
909   *sub_schurs = schurs_ctx;
910   PetscFunctionReturn(0);
911 }
912 
913 #undef __FUNCT__
914 #define __FUNCT__ "PCBDDCSubSchursDestroy"
915 PetscErrorCode PCBDDCSubSchursDestroy(PCBDDCSubSchurs *sub_schurs)
916 {
917   PetscErrorCode ierr;
918 
919   PetscFunctionBegin;
920   ierr = PCBDDCSubSchursReset(*sub_schurs);CHKERRQ(ierr);
921   ierr = PetscFree(*sub_schurs);CHKERRQ(ierr);
922   PetscFunctionReturn(0);
923 }
924 
925 #undef __FUNCT__
926 #define __FUNCT__ "PCBDDCSubSchursReset"
927 PetscErrorCode PCBDDCSubSchursReset(PCBDDCSubSchurs sub_schurs)
928 {
929   PetscInt       i;
930   PetscErrorCode ierr;
931 
932   PetscFunctionBegin;
933   ierr = MatDestroy(&sub_schurs->A);CHKERRQ(ierr);
934   ierr = MatDestroy(&sub_schurs->S);CHKERRQ(ierr);
935   ierr = ISDestroy(&sub_schurs->is_I);CHKERRQ(ierr);
936   ierr = ISDestroy(&sub_schurs->is_B);CHKERRQ(ierr);
937   ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->l2gmap);CHKERRQ(ierr);
938   ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->BtoNmap);CHKERRQ(ierr);
939   ierr = MatDestroy(&sub_schurs->S_Ej_all);CHKERRQ(ierr);
940   ierr = MatDestroy(&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
941   ierr = MatDestroy(&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
942   ierr = MatDestroy(&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
943   ierr = ISDestroy(&sub_schurs->is_Ej_all);CHKERRQ(ierr);
944   ierr = ISDestroy(&sub_schurs->is_Ej_com);CHKERRQ(ierr);
945   ierr = PetscBTDestroy(&sub_schurs->is_edge);CHKERRQ(ierr);
946   ierr = PetscBTDestroy(&sub_schurs->computed_Stilda_subs);CHKERRQ(ierr);
947   for (i=0;i<sub_schurs->n_subs;i++) {
948     ierr = ISDestroy(&sub_schurs->is_subs[i]);CHKERRQ(ierr);
949   }
950   ierr = ISDestroy(&sub_schurs->is_I_layer);CHKERRQ(ierr);
951   if (sub_schurs->n_subs) {
952     ierr = PetscFree(sub_schurs->is_subs);CHKERRQ(ierr);
953   }
954   ierr = PCDestroy(&sub_schurs->interior_solver);CHKERRQ(ierr);
955   sub_schurs->n_subs = 0;
956   PetscFunctionReturn(0);
957 }
958 
959 #undef __FUNCT__
960 #define __FUNCT__ "PCBDDCAdjGetNextLayer_Private"
961 PETSC_STATIC_INLINE PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt* queue_tip,PetscInt n_prev,PetscBT touched,PetscInt* xadj,PetscInt* adjncy,PetscInt* n_added)
962 {
963   PetscInt       i,j,n;
964   PetscErrorCode ierr;
965 
966   PetscFunctionBegin;
967   n = 0;
968   for (i=-n_prev;i<0;i++) {
969     PetscInt start_dof = queue_tip[i];
970     for (j=xadj[start_dof];j<xadj[start_dof+1];j++) {
971       PetscInt dof = adjncy[j];
972       if (!PetscBTLookup(touched,dof)) {
973         ierr = PetscBTSet(touched,dof);CHKERRQ(ierr);
974         queue_tip[n] = dof;
975         n++;
976       }
977     }
978   }
979   *n_added = n;
980   PetscFunctionReturn(0);
981 }
982