xref: /petsc/src/ksp/pc/impls/bddc/bddcschurs.c (revision eb595f79196f12f6b379f375fb4dac46f4df4e83)
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, 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   /* preliminary checks */
211   if (!sub_schurs->use_mumps && compute_Stilda) {
212     SETERRQ(PetscObjectComm((PetscObject)sub_schurs->l2gmap),PETSC_ERR_SUP,"Adaptive selection of constraints requires MUMPS");
213   }
214   /* determine if we are dealing with hermitian positive definite problems */
215   sub_schurs->is_hermitian = PETSC_FALSE;
216   sub_schurs->is_posdef = PETSC_FALSE;
217   if (sub_schurs->A) {
218     PetscInt lsize;
219 
220     ierr = MatGetSize(sub_schurs->A,&lsize,NULL);CHKERRQ(ierr);
221     if (lsize) {
222       ierr = MatIsHermitian(sub_schurs->A,0.0,&sub_schurs->is_hermitian);CHKERRQ(ierr);
223       if (sub_schurs->is_hermitian) {
224         PetscScalar val;
225         Vec         vec1,vec2;
226 
227         ierr = MatCreateVecs(sub_schurs->A,&vec1,&vec2);CHKERRQ(ierr);
228         ierr = VecSetRandom(vec1,NULL);
229         ierr = VecCopy(vec1,vec2);CHKERRQ(ierr);
230         ierr = MatMult(sub_schurs->A,vec2,vec1);CHKERRQ(ierr);
231         ierr = VecDot(vec1,vec2,&val);CHKERRQ(ierr);
232         if (PetscRealPart(val) > 0. && PetscImaginaryPart(val) == 0.) sub_schurs->is_posdef = PETSC_TRUE;
233         ierr = VecDestroy(&vec1);CHKERRQ(ierr);
234         ierr = VecDestroy(&vec2);CHKERRQ(ierr);
235       }
236     } else {
237       sub_schurs->is_hermitian = PETSC_TRUE;
238       sub_schurs->is_posdef = PETSC_TRUE;
239     }
240     if (compute_Stilda && (!sub_schurs->is_hermitian || !sub_schurs->is_posdef)) {
241       SETERRQ(PetscObjectComm((PetscObject)sub_schurs->l2gmap),PETSC_ERR_SUP,"General matrix pencils are not currently supported");
242     }
243   }
244   /* restrict work on active processes */
245   color = 0;
246   if (!sub_schurs->n_subs) color = 1; /* this can happen if we are in a multilevel case or if the subdomain is disconnected */
247   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&rank);CHKERRQ(ierr);
248   ierr = PetscSubcommCreate(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&subcomm);CHKERRQ(ierr);
249   ierr = PetscSubcommSetNumber(subcomm,2);CHKERRQ(ierr);
250   ierr = PetscSubcommSetTypeGeneral(subcomm,color,rank);CHKERRQ(ierr);
251   ierr = PetscCommDuplicate(PetscSubcommChild(subcomm),&comm_n,NULL);CHKERRQ(ierr);
252   ierr = PetscSubcommDestroy(&subcomm);CHKERRQ(ierr);
253   if (!sub_schurs->n_subs) {
254     ierr = PetscCommDestroy(&comm_n);CHKERRQ(ierr);
255     PetscFunctionReturn(0);
256   }
257 
258   /* get Schur complement matrices */
259   if (!sub_schurs->use_mumps) {
260     PetscBool isseqaij;
261     ierr = MatSchurComplementGetSubMatrices(sub_schurs->S,&A_II,NULL,&A_IB,&A_BI,&A_BB);CHKERRQ(ierr);
262     ierr = PetscObjectTypeCompare((PetscObject)A_BB,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
263     if (!isseqaij) {
264       ierr = MatConvert(A_BB,MATSEQAIJ,MAT_REUSE_MATRIX,&A_BB);CHKERRQ(ierr);
265       ierr = MatConvert(A_IB,MATSEQAIJ,MAT_REUSE_MATRIX,&A_IB);CHKERRQ(ierr);
266       ierr = MatConvert(A_BI,MATSEQAIJ,MAT_REUSE_MATRIX,&A_BI);CHKERRQ(ierr);
267     }
268   } else {
269     A_II = NULL;
270     A_IB = NULL;
271     A_BI = NULL;
272     A_BB = NULL;
273   }
274   S_all = NULL;
275   S_all_inv = NULL;
276   S_Ej_tilda_all = NULL;
277   S_Ej_inv_all = NULL;
278 
279   /* determine interior problems */
280   ierr = ISDestroy(&sub_schurs->is_I_layer);CHKERRQ(ierr);
281   ierr = ISGetLocalSize(sub_schurs->is_I,&i);CHKERRQ(ierr);
282   if (nlayers >= 0 && i) { /* Interior problems can be different from the original one */
283     PetscBT                touched;
284     const PetscInt*        idx_B;
285     PetscInt               n_I,n_B,n_local_dofs,n_prev_added,j,layer,*local_numbering;
286 
287     if (xadj == NULL || adjncy == NULL) {
288       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cannot request layering without adjacency");
289     }
290     /* get sizes */
291     ierr = ISGetLocalSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr);
292     ierr = ISGetLocalSize(sub_schurs->is_B,&n_B);CHKERRQ(ierr);
293 
294     ierr = PetscMalloc1(n_I+n_B,&local_numbering);CHKERRQ(ierr);
295     ierr = PetscBTCreate(n_I+n_B,&touched);CHKERRQ(ierr);
296     ierr = PetscBTMemzero(n_I+n_B,touched);CHKERRQ(ierr);
297 
298     /* all boundary dofs must be skipped when adding layers */
299     ierr = ISGetIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr);
300     for (j=0;j<n_B;j++) {
301       ierr = PetscBTSet(touched,idx_B[j]);CHKERRQ(ierr);
302     }
303     ierr = PetscMemcpy(local_numbering,idx_B,n_B*sizeof(PetscInt));CHKERRQ(ierr);
304     ierr = ISRestoreIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr);
305 
306     /* add prescribed number of layers of dofs */
307     n_local_dofs = n_B;
308     n_prev_added = n_B;
309     for (layer=0;layer<nlayers;layer++) {
310       PetscInt n_added;
311       if (n_local_dofs == n_I+n_B) break;
312       if (n_local_dofs > n_I+n_B) {
313         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);
314       }
315       ierr = PCBDDCAdjGetNextLayer_Private(local_numbering+n_local_dofs,n_prev_added,touched,xadj,adjncy,&n_added);CHKERRQ(ierr);
316       n_prev_added = n_added;
317       n_local_dofs += n_added;
318       if (!n_added) break;
319     }
320     ierr = PetscBTDestroy(&touched);CHKERRQ(ierr);
321 
322     /* IS for I layer dofs in original numbering */
323     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);
324     ierr = PetscFree(local_numbering);CHKERRQ(ierr);
325     ierr = ISSort(sub_schurs->is_I_layer);CHKERRQ(ierr);
326     /* IS for I layer dofs in I numbering */
327     if (!sub_schurs->use_mumps) {
328       ISLocalToGlobalMapping ItoNmap;
329       ierr = ISLocalToGlobalMappingCreateIS(sub_schurs->is_I,&ItoNmap);CHKERRQ(ierr);
330       ierr = ISGlobalToLocalMappingApplyIS(ItoNmap,IS_GTOLM_DROP,sub_schurs->is_I_layer,&is_I);CHKERRQ(ierr);
331       ierr = ISLocalToGlobalMappingDestroy(&ItoNmap);CHKERRQ(ierr);
332 
333       /* II block */
334       ierr = MatGetSubMatrix(A_II,is_I,is_I,MAT_INITIAL_MATRIX,&AE_II);CHKERRQ(ierr);
335     }
336   } else {
337     PetscInt n_I;
338 
339     /* IS for I dofs in original numbering */
340     ierr = PetscObjectReference((PetscObject)sub_schurs->is_I);CHKERRQ(ierr);
341     sub_schurs->is_I_layer = sub_schurs->is_I;
342 
343     /* IS for I dofs in I numbering (strided 1) */
344     if (!sub_schurs->use_mumps) {
345       ierr = ISGetSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr);
346       ierr = ISCreateStride(PetscObjectComm((PetscObject)sub_schurs->is_I),n_I,0,1,&is_I);CHKERRQ(ierr);
347 
348       /* II block is the same */
349       ierr = PetscObjectReference((PetscObject)A_II);CHKERRQ(ierr);
350       AE_II = A_II;
351     }
352   }
353 
354   /* Get info on subset sizes and sum of all subsets sizes */
355   max_subset_size = 0;
356   local_size = 0;
357   for (i=0;i<sub_schurs->n_subs;i++) {
358     ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
359     max_subset_size = PetscMax(subset_size,max_subset_size);
360     local_size += subset_size;
361   }
362 
363   /* Work arrays for local indices */
364   extra = 0;
365   if (sub_schurs->use_mumps) {
366     ierr = ISGetLocalSize(sub_schurs->is_I_layer,&extra);CHKERRQ(ierr);
367   }
368   ierr = PetscMalloc1(local_size+extra,&all_local_idx_N);CHKERRQ(ierr);
369   if (extra) {
370     const PetscInt *idxs;
371     ierr = ISGetIndices(sub_schurs->is_I_layer,&idxs);CHKERRQ(ierr);
372     ierr = PetscMemcpy(all_local_idx_N,idxs,extra*sizeof(PetscInt));CHKERRQ(ierr);
373     ierr = ISRestoreIndices(sub_schurs->is_I_layer,&idxs);CHKERRQ(ierr);
374   }
375   ierr = PetscMalloc1(local_size,&nnz);CHKERRQ(ierr);
376   ierr = PetscMalloc2(sub_schurs->n_subs,&auxnum1,sub_schurs->n_subs,&auxnum2);CHKERRQ(ierr);
377 
378   /* Get local indices in local numbering */
379   local_size = 0;
380   for (i=0;i<sub_schurs->n_subs;i++) {
381     PetscInt j;
382     const    PetscInt *idxs;
383 
384     ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
385     ierr = ISGetIndices(sub_schurs->is_subs[i],&idxs);CHKERRQ(ierr);
386     /* start (smallest in global ordering) and multiplicity */
387     auxnum1[i] = idxs[0];
388     auxnum2[i] = subset_size;
389     /* subset indices in local numbering */
390     ierr = PetscMemcpy(all_local_idx_N+local_size+extra,idxs,subset_size*sizeof(PetscInt));CHKERRQ(ierr);
391     ierr = ISRestoreIndices(sub_schurs->is_subs[i],&idxs);CHKERRQ(ierr);
392     for (j=0;j<subset_size;j++) nnz[local_size+j] = subset_size;
393     local_size += subset_size;
394   }
395 
396   /* allocate extra workspace needed only for GETRI */
397   Bwork = NULL;
398   pivots = NULL;
399   if (sub_schurs->n_subs && !sub_schurs->is_hermitian) {
400     PetscScalar lwork;
401 
402     B_lwork = -1;
403     ierr = PetscBLASIntCast(local_size,&B_N);CHKERRQ(ierr);
404     ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
405     PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,Bwork,&B_N,pivots,&lwork,&B_lwork,&B_ierr));
406     ierr = PetscFPTrapPop();CHKERRQ(ierr);
407     if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to GETRI Lapack routine %d",(int)B_ierr);
408     ierr = PetscBLASIntCast((PetscInt)PetscRealPart(lwork),&B_lwork);CHKERRQ(ierr);
409     ierr = PetscMalloc2(B_lwork,&Bwork,B_N,&pivots);CHKERRQ(ierr);
410   }
411 
412   /* prepare parallel matrices for summing up properly schurs on subsets */
413   ierr = PCBDDCSubsetNumbering(comm_n,sub_schurs->l2gmap,sub_schurs->n_subs,auxnum1,auxnum2,&global_size,&all_local_idx_G_rep);CHKERRQ(ierr);
414   ierr = PetscMalloc1(local_size,&all_local_idx_G);CHKERRQ(ierr);
415   local_size = 0;
416   for (i=0;i<sub_schurs->n_subs;i++) {
417     PetscInt j;
418     for (j=0;j<auxnum2[i];j++) all_local_idx_G[local_size++] = all_local_idx_G_rep[i] + j;
419   }
420   ierr = PetscFree(all_local_idx_G_rep);CHKERRQ(ierr);
421   ierr = PetscFree2(auxnum1,auxnum2);CHKERRQ(ierr);
422   ierr = ISLocalToGlobalMappingCreate(comm_n,1,local_size,all_local_idx_G,PETSC_COPY_VALUES,&l2gmap_subsets);CHKERRQ(ierr);
423   ierr = MatCreateIS(comm_n,1,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size,l2gmap_subsets,&work_mat);CHKERRQ(ierr);
424   ierr = ISLocalToGlobalMappingDestroy(&l2gmap_subsets);CHKERRQ(ierr);
425   ierr = MatCreate(PetscObjectComm((PetscObject)work_mat),&global_schur_subsets);CHKERRQ(ierr);
426   ierr = MatSetSizes(global_schur_subsets,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size);CHKERRQ(ierr);
427   ierr = MatSetType(global_schur_subsets,MATMPIAIJ);CHKERRQ(ierr);
428 
429   /* subset indices in local boundary numbering */
430   if (!sub_schurs->is_Ej_all) {
431     PetscInt *all_local_idx_B;
432 
433     ierr = PetscMalloc1(local_size,&all_local_idx_B);CHKERRQ(ierr);
434     ierr = ISGlobalToLocalMappingApply(sub_schurs->BtoNmap,IS_GTOLM_DROP,local_size,all_local_idx_N+extra,&subset_size,all_local_idx_B);CHKERRQ(ierr);
435     if (subset_size != local_size) {
436       SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in sub_schurs serial (BtoNmap)! %d != %d\n",subset_size,local_size);
437     }
438     ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size,all_local_idx_B,PETSC_OWN_POINTER,&sub_schurs->is_Ej_all);CHKERRQ(ierr);
439   }
440 
441   /* Local matrix of all local Schur on subsets (transposed) */
442   if (!sub_schurs->S_Ej_all) {
443     ierr = MatCreate(PETSC_COMM_SELF,&sub_schurs->S_Ej_all);CHKERRQ(ierr);
444     ierr = MatSetSizes(sub_schurs->S_Ej_all,PETSC_DECIDE,PETSC_DECIDE,local_size,local_size);CHKERRQ(ierr);
445     ierr = MatSetType(sub_schurs->S_Ej_all,MATAIJ);CHKERRQ(ierr);
446     ierr = MatSeqAIJSetPreallocation(sub_schurs->S_Ej_all,0,nnz);CHKERRQ(ierr);
447   } else {
448     ierr = MatZeroEntries(sub_schurs->S_Ej_all);CHKERRQ(ierr);
449   }
450 
451   /* Compute Schur complements explicitly */
452   ierr = PetscBTMemzero(sub_schurs->n_subs,sub_schurs->computed_Stilda_subs);CHKERRQ(ierr);
453   if (!sub_schurs->use_mumps) {
454     Mat         S_Ej_expl;
455     PetscScalar *work;
456     PetscInt    j,*dummy_idx;
457     PetscBool   Sdense;
458 
459     ierr = PetscMalloc2(max_subset_size,&dummy_idx,max_subset_size*max_subset_size,&work);CHKERRQ(ierr);
460     local_size = 0;
461     for (i=0;i<sub_schurs->n_subs;i++) {
462       IS  is_subset_B;
463       Mat AE_EE,AE_IE,AE_EI,S_Ej;
464 
465       /* subsets in original and boundary numbering */
466       ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,sub_schurs->is_subs[i],&is_subset_B);CHKERRQ(ierr);
467       /* EE block */
468       ierr = MatGetSubMatrix(A_BB,is_subset_B,is_subset_B,MAT_INITIAL_MATRIX,&AE_EE);CHKERRQ(ierr);
469       /* IE block */
470       ierr = MatGetSubMatrix(A_IB,is_I,is_subset_B,MAT_INITIAL_MATRIX,&AE_IE);CHKERRQ(ierr);
471       /* EI block */
472       if (sub_schurs->is_hermitian) {
473         ierr = MatCreateTranspose(AE_IE,&AE_EI);CHKERRQ(ierr);
474       } else {
475         ierr = MatGetSubMatrix(A_BI,is_subset_B,is_I,MAT_INITIAL_MATRIX,&AE_EI);CHKERRQ(ierr);
476       }
477       ierr = MatCreateSchurComplement(AE_II,AE_II,AE_IE,AE_EI,AE_EE,&S_Ej);CHKERRQ(ierr);
478       ierr = MatDestroy(&AE_EE);CHKERRQ(ierr);
479       ierr = MatDestroy(&AE_IE);CHKERRQ(ierr);
480       ierr = MatDestroy(&AE_EI);CHKERRQ(ierr);
481       if (AE_II == A_II) { /* we can reuse the same ksp */
482         KSP ksp;
483         ierr = MatSchurComplementGetKSP(sub_schurs->S,&ksp);CHKERRQ(ierr);
484         ierr = MatSchurComplementSetKSP(S_Ej,ksp);CHKERRQ(ierr);
485       } else { /* build new ksp object which inherits ksp and pc types from the original one */
486         KSP       origksp,schurksp;
487         PC        origpc,schurpc;
488         KSPType   ksp_type;
489         PetscInt  n_internal;
490         PetscBool ispcnone;
491 
492         ierr = MatSchurComplementGetKSP(sub_schurs->S,&origksp);CHKERRQ(ierr);
493         ierr = MatSchurComplementGetKSP(S_Ej,&schurksp);CHKERRQ(ierr);
494         ierr = KSPGetType(origksp,&ksp_type);CHKERRQ(ierr);
495         ierr = KSPSetType(schurksp,ksp_type);CHKERRQ(ierr);
496         ierr = KSPGetPC(schurksp,&schurpc);CHKERRQ(ierr);
497         ierr = KSPGetPC(origksp,&origpc);CHKERRQ(ierr);
498         ierr = PetscObjectTypeCompare((PetscObject)origpc,PCNONE,&ispcnone);CHKERRQ(ierr);
499         if (!ispcnone) {
500           PCType pc_type;
501           ierr = PCGetType(origpc,&pc_type);CHKERRQ(ierr);
502           ierr = PCSetType(schurpc,pc_type);CHKERRQ(ierr);
503         } else {
504           ierr = PCSetType(schurpc,PCLU);CHKERRQ(ierr);
505         }
506         ierr = ISGetSize(is_I,&n_internal);CHKERRQ(ierr);
507         if (n_internal) { /* UMFPACK gives error with 0 sized problems */
508           MatSolverPackage solver=NULL;
509           ierr = PCFactorGetMatSolverPackage(origpc,(const MatSolverPackage*)&solver);CHKERRQ(ierr);
510           if (solver) {
511             ierr = PCFactorSetMatSolverPackage(schurpc,solver);CHKERRQ(ierr);
512           }
513         }
514         ierr = KSPSetUp(schurksp);CHKERRQ(ierr);
515       }
516       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
517       ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&S_Ej_expl);CHKERRQ(ierr);
518       ierr = PCBDDCComputeExplicitSchur(S_Ej,sub_schurs->is_hermitian,MAT_REUSE_MATRIX,&S_Ej_expl);CHKERRQ(ierr);
519       ierr = PetscObjectTypeCompare((PetscObject)S_Ej_expl,MATSEQDENSE,&Sdense);CHKERRQ(ierr);
520       if (Sdense) {
521         for (j=0;j<subset_size;j++) {
522           dummy_idx[j]=local_size+j;
523         }
524         ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
525       } else {
526         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not yet implemented for sparse matrices");
527       }
528       ierr = MatDestroy(&S_Ej);CHKERRQ(ierr);
529       local_size += subset_size;
530     }
531     ierr = PetscFree2(dummy_idx,work);CHKERRQ(ierr);
532     /* free */
533     ierr = ISDestroy(&is_I);CHKERRQ(ierr);
534     ierr = MatDestroy(&AE_II);CHKERRQ(ierr);
535     ierr = PetscFree(all_local_idx_N);CHKERRQ(ierr);
536   } else {
537     Mat         A,F;
538     IS          is_A_all;
539     PetscScalar *work;
540     PetscInt    *idxs_schur,n_I,n_I_all,*dummy_idx;
541 
542     /* get working mat */
543     ierr = ISGetLocalSize(sub_schurs->is_I_layer,&n_I);CHKERRQ(ierr);
544     ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size+n_I,all_local_idx_N,PETSC_COPY_VALUES,&is_A_all);CHKERRQ(ierr);
545     ierr = MatGetSubMatrixUnsorted(sub_schurs->A,is_A_all,is_A_all,&A);CHKERRQ(ierr);
546     ierr = ISDestroy(&is_A_all);CHKERRQ(ierr);
547 
548     if (n_I) {
549       if (sub_schurs->is_hermitian && sub_schurs->is_posdef) {
550         ierr = MatGetFactor(A,MATSOLVERMUMPS,MAT_FACTOR_CHOLESKY,&F);CHKERRQ(ierr);
551       } else {
552         ierr = MatGetFactor(A,MATSOLVERMUMPS,MAT_FACTOR_LU,&F);CHKERRQ(ierr);
553       }
554 
555       /* subsets ordered last */
556       ierr = PetscMalloc1(local_size,&idxs_schur);CHKERRQ(ierr);
557       for (i=0;i<local_size;i++) {
558         idxs_schur[i] = n_I+i+1;
559       }
560 #if defined(PETSC_HAVE_MUMPS)
561       ierr = MatMumpsSetSchurIndices(F,local_size,idxs_schur);CHKERRQ(ierr);
562 #endif
563       ierr = PetscFree(idxs_schur);CHKERRQ(ierr);
564 
565       /* factorization step */
566       if (sub_schurs->is_hermitian && sub_schurs->is_posdef) {
567         ierr = MatCholeskyFactorSymbolic(F,A,NULL,NULL);CHKERRQ(ierr);
568         ierr = MatCholeskyFactorNumeric(F,A,NULL);CHKERRQ(ierr);
569       } else {
570         ierr = MatLUFactorSymbolic(F,A,NULL,NULL,NULL);CHKERRQ(ierr);
571         ierr = MatLUFactorNumeric(F,A,NULL);CHKERRQ(ierr);
572       }
573 
574       /* get explicit Schur Complement computed during numeric factorization */
575 #if defined(PETSC_HAVE_MUMPS)
576       ierr = MatMumpsGetSchurComplement(F,&S_all);CHKERRQ(ierr);
577 #endif
578 
579       /* we can reuse the interior solver if we are not using the economic version */
580       ierr = ISGetLocalSize(sub_schurs->is_I,&n_I_all);CHKERRQ(ierr);
581       if (n_I == n_I_all) {
582         PCBDDCMumpsInterior msolv_ctx;
583 
584         ierr = PetscNew(&msolv_ctx);CHKERRQ(ierr);
585         msolv_ctx->n = n_I;
586         ierr = PetscObjectReference((PetscObject)F);CHKERRQ(ierr);
587         msolv_ctx->F = F;
588         ierr = MatCreateVecs(F,&msolv_ctx->sol,&msolv_ctx->rhs);CHKERRQ(ierr);
589         ierr = PCCreate(PETSC_COMM_SELF,&sub_schurs->interior_solver);CHKERRQ(ierr);
590         ierr = MatSchurComplementGetSubMatrices(sub_schurs->S,&A_II,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
591         ierr = PCSetOperators(sub_schurs->interior_solver,A_II,A_II);CHKERRQ(ierr);
592         ierr = PCSetType(sub_schurs->interior_solver,PCSHELL);CHKERRQ(ierr);
593         ierr = PCShellSetContext(sub_schurs->interior_solver,msolv_ctx);CHKERRQ(ierr);
594         ierr = PCShellSetApply(sub_schurs->interior_solver,PCBDDCMumpsInteriorSolve);CHKERRQ(ierr);
595         ierr = PCShellSetDestroy(sub_schurs->interior_solver,PCBDDCMumpsInteriorDestroy);CHKERRQ(ierr);
596       }
597       ierr = MatDestroy(&F);CHKERRQ(ierr);
598     } else {
599       ierr = MatConvert(A,MATSEQDENSE,MAT_INITIAL_MATRIX,&S_all);CHKERRQ(ierr);
600     }
601     ierr = MatDestroy(&A);CHKERRQ(ierr);
602     ierr = PetscFree(all_local_idx_N);CHKERRQ(ierr);
603 
604     if (compute_Stilda) {
605       ierr = MatCreate(PETSC_COMM_SELF,&S_Ej_tilda_all);CHKERRQ(ierr);
606       ierr = MatSetSizes(S_Ej_tilda_all,PETSC_DECIDE,PETSC_DECIDE,local_size,local_size);CHKERRQ(ierr);
607       ierr = MatSetType(S_Ej_tilda_all,MATAIJ);CHKERRQ(ierr);
608       ierr = MatSeqAIJSetPreallocation(S_Ej_tilda_all,0,nnz);CHKERRQ(ierr);
609       ierr = MatCreate(PETSC_COMM_SELF,&S_Ej_inv_all);CHKERRQ(ierr);
610       ierr = MatSetSizes(S_Ej_inv_all,PETSC_DECIDE,PETSC_DECIDE,local_size,local_size);CHKERRQ(ierr);
611       ierr = MatSetType(S_Ej_inv_all,MATAIJ);CHKERRQ(ierr);
612       ierr = MatSeqAIJSetPreallocation(S_Ej_inv_all,0,nnz);CHKERRQ(ierr);
613 
614       /* compute St^-1 */
615       if (local_size) { /* multilevel guard */
616         PetscScalar *vals;
617 
618         ierr = PetscBLASIntCast(local_size,&B_N);CHKERRQ(ierr);
619         ierr = MatDuplicate(S_all,MAT_COPY_VALUES,&S_all_inv);CHKERRQ(ierr);
620         ierr = MatDenseGetArray(S_all_inv,&vals);CHKERRQ(ierr);
621         if (!sub_schurs->is_hermitian) {
622           PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,vals,&B_N,pivots,&B_ierr));
623           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
624           PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,vals,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
625           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
626         } else {
627           PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,vals,&B_N,&B_ierr));
628           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
629           PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,vals,&B_N,&B_ierr));
630           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
631         }
632         ierr = MatDenseRestoreArray(S_all_inv,&vals);CHKERRQ(ierr);
633       }
634     }
635 
636     /* Work arrays */
637     if (sub_schurs->n_subs == 1) {
638       ierr = PetscMalloc1(max_subset_size,&dummy_idx);CHKERRQ(ierr);
639     } else {
640       ierr = PetscMalloc2(max_subset_size,&dummy_idx,max_subset_size*max_subset_size,&work);CHKERRQ(ierr);
641     }
642 
643     local_size = 0;
644     for (i=0;i<sub_schurs->n_subs;i++) {
645       Mat S_Ej;
646       IS  is_E;
647       PetscInt j;
648 
649       /* get S_E */
650       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
651       if (sub_schurs->n_subs == 1) {
652         ierr = MatDenseGetArray(S_all,&work);CHKERRQ(ierr);
653         S_Ej = NULL;
654         is_E = NULL;
655       } else {
656         ierr = ISCreateStride(PETSC_COMM_SELF,subset_size,local_size,1,&is_E);CHKERRQ(ierr);
657         ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&S_Ej);CHKERRQ(ierr);
658         ierr = MatGetSubMatrix(S_all,is_E,is_E,MAT_REUSE_MATRIX,&S_Ej);CHKERRQ(ierr);
659       }
660       /* insert S_E values */
661       for (j=0;j<subset_size;j++) {
662         dummy_idx[j]=local_size+j;
663       }
664       ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
665 
666       /* if adaptivity is requested, invert S_E and insert St_E^-1 blocks */
667       if (compute_Stilda && ((PetscBTLookup(sub_schurs->is_edge,i) && use_edges) || (!PetscBTLookup(sub_schurs->is_edge,i) && use_faces))) {
668         /* get S_E^-1 */
669         ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr);
670         ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
671         if (!sub_schurs->is_hermitian) {
672           PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,work,&B_N,pivots,&B_ierr));
673           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
674           PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,work,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
675           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
676         } else {
677           PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,work,&B_N,&B_ierr));
678           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
679           PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,work,&B_N,&B_ierr));
680           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
681         }
682         ierr = PetscFPTrapPop();CHKERRQ(ierr);
683         ierr = MatSetValues(S_Ej_inv_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
684 
685         /* get St_E^-1 */
686         if (sub_schurs->n_subs == 1) {
687           ierr = MatDenseRestoreArray(S_all,&work);CHKERRQ(ierr);
688           ierr = MatDenseGetArray(S_all_inv,&work);CHKERRQ(ierr);
689         } else {
690           ierr = MatGetSubMatrix(S_all_inv,is_E,is_E,MAT_REUSE_MATRIX,&S_Ej);CHKERRQ(ierr);
691         }
692         ierr = MatSetValues(S_Ej_tilda_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
693         if (sub_schurs->n_subs == 1) {
694           ierr = MatDenseRestoreArray(S_all_inv,&work);CHKERRQ(ierr);
695         }
696         ierr = PetscBTSet(sub_schurs->computed_Stilda_subs,i);CHKERRQ(ierr);
697       } else if (sub_schurs->n_subs == 1) {
698         ierr = MatDenseRestoreArray(S_all,&work);CHKERRQ(ierr);
699       }
700       ierr = MatDestroy(&S_Ej);CHKERRQ(ierr);
701       ierr = ISDestroy(&is_E);CHKERRQ(ierr);
702       local_size += subset_size;
703     }
704     if (sub_schurs->n_subs == 1) {
705       ierr = PetscFree(dummy_idx);CHKERRQ(ierr);
706     } else {
707       ierr = PetscFree2(dummy_idx,work);CHKERRQ(ierr);
708     }
709   }
710   ierr = PetscFree(nnz);CHKERRQ(ierr);
711   ierr = MatDestroy(&S_all);CHKERRQ(ierr);
712   ierr = MatDestroy(&S_all_inv);CHKERRQ(ierr);
713   ierr = MatAssemblyBegin(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
714   ierr = MatAssemblyEnd(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
715   if (compute_Stilda) {
716     ierr = MatAssemblyBegin(S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
717     ierr = MatAssemblyEnd(S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
718     ierr = MatAssemblyBegin(S_Ej_inv_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
719     ierr = MatAssemblyEnd(S_Ej_inv_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
720   }
721 
722   /* Global matrix of all assembled Schur on subsets */
723   ierr = MatISSetLocalMat(work_mat,sub_schurs->S_Ej_all);CHKERRQ(ierr);
724   ierr = MatISSetMPIXAIJPreallocation_Private(work_mat,global_schur_subsets,PETSC_TRUE);CHKERRQ(ierr);
725   ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
726 
727   /* Get local part of (\sum_j S_Ej) */
728   ierr = ISCreateGeneral(comm_n,local_size,all_local_idx_G,PETSC_OWN_POINTER,&temp_is);CHKERRQ(ierr);
729   ierr = MatDestroy(&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
730   ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
731 
732   /* Compute explicitly (\sum_j S_Ej)^-1 (faster scaling during PCApply, needs extra work when doing setup) */
733   if (faster_deluxe) {
734     Mat         tmpmat;
735     PetscScalar *array;
736     PetscInt    cum;
737 
738     ierr = MatSeqAIJGetArray(sub_schurs->sum_S_Ej_all,&array);CHKERRQ(ierr);
739     cum = 0;
740     for (i=0;i<sub_schurs->n_subs;i++) {
741       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
742       ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr);
743       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
744       if (!sub_schurs->is_hermitian) {
745         PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,array+cum,&B_N,pivots,&B_ierr));
746         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
747         PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,array+cum,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
748         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
749       } else {
750         PetscInt j,k;
751 
752         PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,array+cum,&B_N,&B_ierr));
753         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
754         PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,array+cum,&B_N,&B_ierr));
755         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
756         for (j=0;j<B_N;j++) {
757           for (k=j+1;k<B_N;k++) {
758             array[k*B_N+j+cum] = array[j*B_N+k+cum];
759           }
760         }
761       }
762       ierr = PetscFPTrapPop();CHKERRQ(ierr);
763       cum += subset_size*subset_size;
764     }
765     ierr = MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_all,&array);CHKERRQ(ierr);
766     ierr = MatMatMult(sub_schurs->S_Ej_all,sub_schurs->sum_S_Ej_all,MAT_INITIAL_MATRIX,1.0,&tmpmat);CHKERRQ(ierr);
767     ierr = MatDestroy(&sub_schurs->S_Ej_all);CHKERRQ(ierr);
768     ierr = MatDestroy(&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
769     sub_schurs->S_Ej_all = tmpmat;
770   }
771 
772   /* Get local part of (\sum_j S^-1_Ej) (\sum_j St^-1_Ej) */
773   if (compute_Stilda) {
774     ierr = MatISSetLocalMat(work_mat,S_Ej_tilda_all);CHKERRQ(ierr);
775     ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
776     ierr = MatDestroy(&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
777     ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
778     ierr = MatISSetLocalMat(work_mat,S_Ej_inv_all);CHKERRQ(ierr);
779     ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
780     ierr = MatDestroy(&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
781     ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
782   }
783 
784   /* free workspace */
785   ierr = PetscFree2(Bwork,pivots);CHKERRQ(ierr);
786   ierr = MatDestroy(&global_schur_subsets);CHKERRQ(ierr);
787   ierr = MatDestroy(&S_Ej_tilda_all);CHKERRQ(ierr);
788   ierr = MatDestroy(&S_Ej_inv_all);CHKERRQ(ierr);
789   ierr = MatDestroy(&work_mat);CHKERRQ(ierr);
790   ierr = ISDestroy(&temp_is);CHKERRQ(ierr);
791   ierr = PetscCommDestroy(&comm_n);CHKERRQ(ierr);
792   PetscFunctionReturn(0);
793 }
794 
795 #undef __FUNCT__
796 #define __FUNCT__ "PCBDDCSubSchursInit"
797 PetscErrorCode PCBDDCSubSchursInit(PCBDDCSubSchurs sub_schurs, Mat A, Mat S, IS is_I, IS is_B, PCBDDCGraph graph, ISLocalToGlobalMapping BtoNmap)
798 {
799   IS              *faces,*edges,*all_cc,vertices;
800   PetscInt        i,n_faces,n_edges,n_all_cc;
801   PetscBool       is_sorted;
802   PetscErrorCode  ierr;
803 
804   PetscFunctionBegin;
805   ierr = ISSorted(is_I,&is_sorted);CHKERRQ(ierr);
806   if (!is_sorted) {
807     SETERRQ(PetscObjectComm((PetscObject)is_I),PETSC_ERR_PLIB,"IS for I dofs should be shorted");
808   }
809   ierr = ISSorted(is_B,&is_sorted);CHKERRQ(ierr);
810   if (!is_sorted) {
811     SETERRQ(PetscObjectComm((PetscObject)is_B),PETSC_ERR_PLIB,"IS for B dofs should be shorted");
812   }
813 
814   /* reset any previous data */
815   ierr = PCBDDCSubSchursReset(sub_schurs);CHKERRQ(ierr);
816 
817   /* get index sets for faces and edges (already sorted by global ordering) */
818   ierr = PCBDDCGraphGetCandidatesIS(graph,&n_faces,&faces,&n_edges,&edges,&vertices);CHKERRQ(ierr);
819   n_all_cc = n_faces+n_edges;
820   ierr = PetscBTCreate(n_all_cc,&sub_schurs->is_edge);CHKERRQ(ierr);
821   ierr = PetscBTCreate(n_all_cc,&sub_schurs->computed_Stilda_subs);CHKERRQ(ierr);
822   ierr = PetscMalloc1(n_all_cc,&all_cc);CHKERRQ(ierr);
823   for (i=0;i<n_faces;i++) {
824     all_cc[i] = faces[i];
825   }
826   for (i=0;i<n_edges;i++) {
827     all_cc[n_faces+i] = edges[i];
828     ierr = PetscBTSet(sub_schurs->is_edge,n_faces+i);CHKERRQ(ierr);
829   }
830   ierr = PetscFree(faces);CHKERRQ(ierr);
831   ierr = PetscFree(edges);CHKERRQ(ierr);
832 
833   /* Determine if MUMPS cane be used */
834   sub_schurs->use_mumps = PETSC_FALSE;
835 #if defined(PETSC_HAVE_MUMPS)
836   sub_schurs->use_mumps = (PetscBool)(!!A);
837 #endif
838 
839   /* update info in sub_schurs */
840   if (A) {
841     PetscBool isseqaij;
842 
843     ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
844     if (isseqaij) {
845       ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
846       sub_schurs->A = A;
847     } else { /* SeqBAIJ matrices does not support symmetry checking, SeqSBAIJ does not support MatPermute */
848       ierr = MatConvert(A,MATSEQAIJ,MAT_INITIAL_MATRIX,&sub_schurs->A);CHKERRQ(ierr);
849     }
850   }
851   ierr = PetscObjectReference((PetscObject)S);CHKERRQ(ierr);
852   sub_schurs->S = S;
853   ierr = PetscObjectReference((PetscObject)is_I);CHKERRQ(ierr);
854   sub_schurs->is_I = is_I;
855   ierr = PetscObjectReference((PetscObject)is_B);CHKERRQ(ierr);
856   sub_schurs->is_B = is_B;
857   ierr = PetscObjectReference((PetscObject)graph->l2gmap);CHKERRQ(ierr);
858   sub_schurs->l2gmap = graph->l2gmap;
859   ierr = PetscObjectReference((PetscObject)BtoNmap);CHKERRQ(ierr);
860   sub_schurs->BtoNmap = BtoNmap;
861   sub_schurs->n_subs = n_all_cc;
862   sub_schurs->is_subs = all_cc;
863   if (!sub_schurs->use_mumps) { /* for adaptive selection */
864     for (i=0;i<sub_schurs->n_subs;i++) {
865       ierr = ISSort(sub_schurs->is_subs[i]);CHKERRQ(ierr);
866     }
867   }
868   sub_schurs->is_Ej_com = vertices;
869 
870 
871   /* allocate space for schur complements */
872   sub_schurs->S_Ej_all = NULL;
873   sub_schurs->sum_S_Ej_all = NULL;
874   sub_schurs->sum_S_Ej_inv_all = NULL;
875   sub_schurs->sum_S_Ej_tilda_all = NULL;
876   sub_schurs->is_Ej_all = NULL;
877   PetscFunctionReturn(0);
878 }
879 
880 #undef __FUNCT__
881 #define __FUNCT__ "PCBDDCSubSchursCreate"
882 PetscErrorCode PCBDDCSubSchursCreate(PCBDDCSubSchurs *sub_schurs)
883 {
884   PCBDDCSubSchurs schurs_ctx;
885   PetscErrorCode  ierr;
886 
887   PetscFunctionBegin;
888   ierr = PetscNew(&schurs_ctx);CHKERRQ(ierr);
889   schurs_ctx->n_subs = 0;
890   *sub_schurs = schurs_ctx;
891   PetscFunctionReturn(0);
892 }
893 
894 #undef __FUNCT__
895 #define __FUNCT__ "PCBDDCSubSchursDestroy"
896 PetscErrorCode PCBDDCSubSchursDestroy(PCBDDCSubSchurs *sub_schurs)
897 {
898   PetscErrorCode ierr;
899 
900   PetscFunctionBegin;
901   ierr = PCBDDCSubSchursReset(*sub_schurs);CHKERRQ(ierr);
902   ierr = PetscFree(*sub_schurs);CHKERRQ(ierr);
903   PetscFunctionReturn(0);
904 }
905 
906 #undef __FUNCT__
907 #define __FUNCT__ "PCBDDCSubSchursReset"
908 PetscErrorCode PCBDDCSubSchursReset(PCBDDCSubSchurs sub_schurs)
909 {
910   PetscInt       i;
911   PetscErrorCode ierr;
912 
913   PetscFunctionBegin;
914   ierr = MatDestroy(&sub_schurs->A);CHKERRQ(ierr);
915   ierr = MatDestroy(&sub_schurs->S);CHKERRQ(ierr);
916   ierr = ISDestroy(&sub_schurs->is_I);CHKERRQ(ierr);
917   ierr = ISDestroy(&sub_schurs->is_B);CHKERRQ(ierr);
918   ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->l2gmap);CHKERRQ(ierr);
919   ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->BtoNmap);CHKERRQ(ierr);
920   ierr = MatDestroy(&sub_schurs->S_Ej_all);CHKERRQ(ierr);
921   ierr = MatDestroy(&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
922   ierr = MatDestroy(&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
923   ierr = MatDestroy(&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
924   ierr = ISDestroy(&sub_schurs->is_Ej_all);CHKERRQ(ierr);
925   ierr = ISDestroy(&sub_schurs->is_Ej_com);CHKERRQ(ierr);
926   ierr = PetscBTDestroy(&sub_schurs->is_edge);CHKERRQ(ierr);
927   ierr = PetscBTDestroy(&sub_schurs->computed_Stilda_subs);CHKERRQ(ierr);
928   for (i=0;i<sub_schurs->n_subs;i++) {
929     ierr = ISDestroy(&sub_schurs->is_subs[i]);CHKERRQ(ierr);
930   }
931   ierr = ISDestroy(&sub_schurs->is_I_layer);CHKERRQ(ierr);
932   if (sub_schurs->n_subs) {
933     ierr = PetscFree(sub_schurs->is_subs);CHKERRQ(ierr);
934   }
935   ierr = PCDestroy(&sub_schurs->interior_solver);CHKERRQ(ierr);
936   sub_schurs->n_subs = 0;
937   PetscFunctionReturn(0);
938 }
939 
940 #undef __FUNCT__
941 #define __FUNCT__ "PCBDDCAdjGetNextLayer_Private"
942 PETSC_STATIC_INLINE PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt* queue_tip,PetscInt n_prev,PetscBT touched,PetscInt* xadj,PetscInt* adjncy,PetscInt* n_added)
943 {
944   PetscInt       i,j,n;
945   PetscErrorCode ierr;
946 
947   PetscFunctionBegin;
948   n = 0;
949   for (i=-n_prev;i<0;i++) {
950     PetscInt start_dof = queue_tip[i];
951     for (j=xadj[start_dof];j<xadj[start_dof+1];j++) {
952       PetscInt dof = adjncy[j];
953       if (!PetscBTLookup(touched,dof)) {
954         ierr = PetscBTSet(touched,dof);CHKERRQ(ierr);
955         queue_tip[n] = dof;
956         n++;
957       }
958     }
959   }
960   *n_added = n;
961   PetscFunctionReturn(0);
962 }
963