xref: /petsc/src/ksp/pc/impls/bddc/bddcschurs.c (revision a9b995523724ee1ed9940043d5e29baebc65db58)
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,is_I_layer,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 = ISGetLocalSize(sub_schurs->is_I,&i);CHKERRQ(ierr);
309   if (nlayers >= 0 && i) { /* Interior problems can be different from the original one */
310     PetscBT                touched;
311     const PetscInt*        idx_B;
312     PetscInt               n_I,n_B,n_local_dofs,n_prev_added,j,layer,*local_numbering;
313 
314     if (xadj == NULL || adjncy == NULL) {
315       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cannot request layering without adjacency");
316     }
317     /* get sizes */
318     ierr = ISGetLocalSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr);
319     ierr = ISGetLocalSize(sub_schurs->is_B,&n_B);CHKERRQ(ierr);
320 
321     ierr = PetscMalloc1(n_I+n_B,&local_numbering);CHKERRQ(ierr);
322     ierr = PetscBTCreate(n_I+n_B,&touched);CHKERRQ(ierr);
323     ierr = PetscBTMemzero(n_I+n_B,touched);CHKERRQ(ierr);
324 
325     /* all boundary dofs must be skipped when adding layers */
326     ierr = ISGetIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr);
327     for (j=0;j<n_B;j++) {
328       ierr = PetscBTSet(touched,idx_B[j]);CHKERRQ(ierr);
329     }
330     ierr = PetscMemcpy(local_numbering,idx_B,n_B*sizeof(PetscInt));CHKERRQ(ierr);
331     ierr = ISRestoreIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr);
332 
333     /* add prescribed number of layers of dofs */
334     n_local_dofs = n_B;
335     n_prev_added = n_B;
336     for (layer=0;layer<nlayers;layer++) {
337       PetscInt n_added;
338       if (n_local_dofs == n_I+n_B) break;
339       if (n_local_dofs > n_I+n_B) {
340         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);
341       }
342       ierr = PCBDDCAdjGetNextLayer_Private(local_numbering+n_local_dofs,n_prev_added,touched,xadj,adjncy,&n_added);CHKERRQ(ierr);
343       n_prev_added = n_added;
344       n_local_dofs += n_added;
345       if (!n_added) break;
346     }
347     ierr = PetscBTDestroy(&touched);CHKERRQ(ierr);
348 
349     /* IS for I layer dofs in original numbering */
350     ierr = ISCreateGeneral(PetscObjectComm((PetscObject)sub_schurs->is_I),n_local_dofs-n_B,local_numbering+n_B,PETSC_COPY_VALUES,&is_I_layer);CHKERRQ(ierr);
351     ierr = PetscFree(local_numbering);CHKERRQ(ierr);
352     ierr = ISSort(is_I_layer);CHKERRQ(ierr);
353     /* IS for I layer dofs in I numbering */
354     if (!sub_schurs->use_mumps) {
355       ISLocalToGlobalMapping ItoNmap;
356       ierr = ISLocalToGlobalMappingCreateIS(sub_schurs->is_I,&ItoNmap);CHKERRQ(ierr);
357       ierr = ISGlobalToLocalMappingApplyIS(ItoNmap,IS_GTOLM_DROP,is_I_layer,&is_I);CHKERRQ(ierr);
358       ierr = ISLocalToGlobalMappingDestroy(&ItoNmap);CHKERRQ(ierr);
359 
360       /* II block */
361       ierr = MatGetSubMatrix(A_II,is_I,is_I,MAT_INITIAL_MATRIX,&AE_II);CHKERRQ(ierr);
362     }
363   } else {
364     PetscInt n_I;
365 
366     /* IS for I dofs in original numbering */
367     ierr = PetscObjectReference((PetscObject)sub_schurs->is_I);CHKERRQ(ierr);
368     is_I_layer = sub_schurs->is_I;
369 
370     /* IS for I dofs in I numbering (strided 1) */
371     if (!sub_schurs->use_mumps) {
372       ierr = ISGetSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr);
373       ierr = ISCreateStride(PetscObjectComm((PetscObject)sub_schurs->is_I),n_I,0,1,&is_I);CHKERRQ(ierr);
374 
375       /* II block is the same */
376       ierr = PetscObjectReference((PetscObject)A_II);CHKERRQ(ierr);
377       AE_II = A_II;
378     }
379   }
380 
381   /* Get info on subset sizes and sum of all subsets sizes */
382   max_subset_size = 0;
383   local_size = 0;
384   for (i=0;i<sub_schurs->n_subs;i++) {
385     ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
386     max_subset_size = PetscMax(subset_size,max_subset_size);
387     local_size += subset_size;
388   }
389 
390   /* Work arrays for local indices */
391   extra = 0;
392   if (sub_schurs->use_mumps) {
393     ierr = ISGetLocalSize(is_I_layer,&extra);CHKERRQ(ierr);
394   }
395   ierr = PetscMalloc1(local_size+extra,&all_local_idx_N);CHKERRQ(ierr);
396   if (extra) {
397     const PetscInt *idxs;
398     ierr = ISGetIndices(is_I_layer,&idxs);CHKERRQ(ierr);
399     ierr = PetscMemcpy(all_local_idx_N,idxs,extra*sizeof(PetscInt));CHKERRQ(ierr);
400     ierr = ISRestoreIndices(is_I_layer,&idxs);CHKERRQ(ierr);
401   }
402   ierr = PetscMalloc1(local_size,&nnz);CHKERRQ(ierr);
403   ierr = PetscMalloc2(sub_schurs->n_subs,&auxnum1,sub_schurs->n_subs,&auxnum2);CHKERRQ(ierr);
404 
405   /* Get local indices in local numbering */
406   local_size = 0;
407   for (i=0;i<sub_schurs->n_subs;i++) {
408     PetscInt j;
409     const    PetscInt *idxs;
410 
411     ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
412     ierr = ISGetIndices(sub_schurs->is_subs[i],&idxs);CHKERRQ(ierr);
413     /* start (smallest in global ordering) and multiplicity */
414     auxnum1[i] = idxs[0];
415     auxnum2[i] = subset_size;
416     /* subset indices in local numbering */
417     ierr = PetscMemcpy(all_local_idx_N+local_size+extra,idxs,subset_size*sizeof(PetscInt));CHKERRQ(ierr);
418     ierr = ISRestoreIndices(sub_schurs->is_subs[i],&idxs);CHKERRQ(ierr);
419     for (j=0;j<subset_size;j++) nnz[local_size+j] = subset_size;
420     local_size += subset_size;
421   }
422 
423   /* allocate extra workspace needed only for GETRI */
424   Bwork = NULL;
425   pivots = NULL;
426   if (sub_schurs->n_subs && !sub_schurs->is_hermitian) {
427     PetscScalar lwork;
428 
429     B_lwork = -1;
430     ierr = PetscBLASIntCast(local_size,&B_N);CHKERRQ(ierr);
431     ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
432     PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,Bwork,&B_N,pivots,&lwork,&B_lwork,&B_ierr));
433     ierr = PetscFPTrapPop();CHKERRQ(ierr);
434     if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to GETRI Lapack routine %d",(int)B_ierr);
435     ierr = PetscBLASIntCast((PetscInt)PetscRealPart(lwork),&B_lwork);CHKERRQ(ierr);
436     ierr = PetscMalloc2(B_lwork,&Bwork,B_N,&pivots);CHKERRQ(ierr);
437   }
438 
439   /* prepare parallel matrices for summing up properly schurs on subsets */
440   ierr = PCBDDCSubsetNumbering(comm_n,sub_schurs->l2gmap,sub_schurs->n_subs,auxnum1,auxnum2,&global_size,&all_local_idx_G_rep);CHKERRQ(ierr);
441   ierr = PetscMalloc1(local_size,&all_local_idx_G);CHKERRQ(ierr);
442   local_size = 0;
443   for (i=0;i<sub_schurs->n_subs;i++) {
444     PetscInt j;
445     for (j=0;j<auxnum2[i];j++) all_local_idx_G[local_size++] = all_local_idx_G_rep[i] + j;
446   }
447   ierr = PetscFree(all_local_idx_G_rep);CHKERRQ(ierr);
448   ierr = PetscFree2(auxnum1,auxnum2);CHKERRQ(ierr);
449   ierr = ISLocalToGlobalMappingCreate(comm_n,1,local_size,all_local_idx_G,PETSC_COPY_VALUES,&l2gmap_subsets);CHKERRQ(ierr);
450   ierr = MatCreateIS(comm_n,1,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size,l2gmap_subsets,&work_mat);CHKERRQ(ierr);
451   ierr = ISLocalToGlobalMappingDestroy(&l2gmap_subsets);CHKERRQ(ierr);
452   ierr = MatCreate(PetscObjectComm((PetscObject)work_mat),&global_schur_subsets);CHKERRQ(ierr);
453   ierr = MatSetSizes(global_schur_subsets,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size);CHKERRQ(ierr);
454   ierr = MatSetType(global_schur_subsets,MATMPIAIJ);CHKERRQ(ierr);
455 
456   /* subset indices in local boundary numbering */
457   if (!sub_schurs->is_Ej_all) {
458     PetscInt *all_local_idx_B;
459 
460     ierr = PetscMalloc1(local_size,&all_local_idx_B);CHKERRQ(ierr);
461     ierr = ISGlobalToLocalMappingApply(sub_schurs->BtoNmap,IS_GTOLM_DROP,local_size,all_local_idx_N+extra,&subset_size,all_local_idx_B);CHKERRQ(ierr);
462     if (subset_size != local_size) {
463       SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in sub_schurs serial (BtoNmap)! %d != %d\n",subset_size,local_size);
464     }
465     ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size,all_local_idx_B,PETSC_OWN_POINTER,&sub_schurs->is_Ej_all);CHKERRQ(ierr);
466   }
467 
468   /* Local matrix of all local Schur on subsets (transposed) */
469   if (!sub_schurs->S_Ej_all) {
470     ierr = MatCreate(PETSC_COMM_SELF,&sub_schurs->S_Ej_all);CHKERRQ(ierr);
471     ierr = MatSetSizes(sub_schurs->S_Ej_all,PETSC_DECIDE,PETSC_DECIDE,local_size,local_size);CHKERRQ(ierr);
472     ierr = MatSetType(sub_schurs->S_Ej_all,MATAIJ);CHKERRQ(ierr);
473     ierr = MatSeqAIJSetPreallocation(sub_schurs->S_Ej_all,0,nnz);CHKERRQ(ierr);
474   } else {
475     ierr = MatZeroEntries(sub_schurs->S_Ej_all);CHKERRQ(ierr);
476   }
477 
478   /* Compute Schur complements explicitly */
479   ierr = PetscBTMemzero(sub_schurs->n_subs,sub_schurs->computed_Stilda_subs);CHKERRQ(ierr);
480   if (!sub_schurs->use_mumps) {
481     Mat         S_Ej_expl;
482     PetscScalar *work;
483     PetscInt    j,*dummy_idx;
484     PetscBool   Sdense;
485 
486     ierr = PetscMalloc2(max_subset_size,&dummy_idx,max_subset_size*max_subset_size,&work);CHKERRQ(ierr);
487     local_size = 0;
488     for (i=0;i<sub_schurs->n_subs;i++) {
489       IS  is_subset_B;
490       Mat AE_EE,AE_IE,AE_EI,S_Ej;
491 
492       /* subsets in original and boundary numbering */
493       ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,sub_schurs->is_subs[i],&is_subset_B);CHKERRQ(ierr);
494       /* EE block */
495       ierr = MatGetSubMatrix(A_BB,is_subset_B,is_subset_B,MAT_INITIAL_MATRIX,&AE_EE);CHKERRQ(ierr);
496       /* IE block */
497       ierr = MatGetSubMatrix(A_IB,is_I,is_subset_B,MAT_INITIAL_MATRIX,&AE_IE);CHKERRQ(ierr);
498       /* EI block */
499       if (sub_schurs->is_hermitian) {
500         ierr = MatCreateTranspose(AE_IE,&AE_EI);CHKERRQ(ierr);
501       } else {
502         ierr = MatGetSubMatrix(A_BI,is_subset_B,is_I,MAT_INITIAL_MATRIX,&AE_EI);CHKERRQ(ierr);
503       }
504       ierr = ISDestroy(&is_subset_B);CHKERRQ(ierr);
505       ierr = MatCreateSchurComplement(AE_II,AE_II,AE_IE,AE_EI,AE_EE,&S_Ej);CHKERRQ(ierr);
506       ierr = MatDestroy(&AE_EE);CHKERRQ(ierr);
507       ierr = MatDestroy(&AE_IE);CHKERRQ(ierr);
508       ierr = MatDestroy(&AE_EI);CHKERRQ(ierr);
509       if (AE_II == A_II) { /* we can reuse the same ksp */
510         KSP ksp;
511         ierr = MatSchurComplementGetKSP(sub_schurs->S,&ksp);CHKERRQ(ierr);
512         ierr = MatSchurComplementSetKSP(S_Ej,ksp);CHKERRQ(ierr);
513       } else { /* build new ksp object which inherits ksp and pc types from the original one */
514         KSP       origksp,schurksp;
515         PC        origpc,schurpc;
516         KSPType   ksp_type;
517         PetscInt  n_internal;
518         PetscBool ispcnone;
519 
520         ierr = MatSchurComplementGetKSP(sub_schurs->S,&origksp);CHKERRQ(ierr);
521         ierr = MatSchurComplementGetKSP(S_Ej,&schurksp);CHKERRQ(ierr);
522         ierr = KSPGetType(origksp,&ksp_type);CHKERRQ(ierr);
523         ierr = KSPSetType(schurksp,ksp_type);CHKERRQ(ierr);
524         ierr = KSPGetPC(schurksp,&schurpc);CHKERRQ(ierr);
525         ierr = KSPGetPC(origksp,&origpc);CHKERRQ(ierr);
526         ierr = PetscObjectTypeCompare((PetscObject)origpc,PCNONE,&ispcnone);CHKERRQ(ierr);
527         if (!ispcnone) {
528           PCType pc_type;
529           ierr = PCGetType(origpc,&pc_type);CHKERRQ(ierr);
530           ierr = PCSetType(schurpc,pc_type);CHKERRQ(ierr);
531         } else {
532           ierr = PCSetType(schurpc,PCLU);CHKERRQ(ierr);
533         }
534         ierr = ISGetSize(is_I,&n_internal);CHKERRQ(ierr);
535         if (n_internal) { /* UMFPACK gives error with 0 sized problems */
536           MatSolverPackage solver=NULL;
537           ierr = PCFactorGetMatSolverPackage(origpc,(const MatSolverPackage*)&solver);CHKERRQ(ierr);
538           if (solver) {
539             ierr = PCFactorSetMatSolverPackage(schurpc,solver);CHKERRQ(ierr);
540           }
541         }
542         ierr = KSPSetUp(schurksp);CHKERRQ(ierr);
543       }
544       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
545       ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&S_Ej_expl);CHKERRQ(ierr);
546       ierr = PCBDDCComputeExplicitSchur(S_Ej,sub_schurs->is_hermitian,MAT_REUSE_MATRIX,&S_Ej_expl);CHKERRQ(ierr);
547       ierr = PetscObjectTypeCompare((PetscObject)S_Ej_expl,MATSEQDENSE,&Sdense);CHKERRQ(ierr);
548       if (Sdense) {
549         for (j=0;j<subset_size;j++) {
550           dummy_idx[j]=local_size+j;
551         }
552         ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
553       } else {
554         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not yet implemented for sparse matrices");
555       }
556       ierr = MatDestroy(&S_Ej);CHKERRQ(ierr);
557       ierr = MatDestroy(&S_Ej_expl);CHKERRQ(ierr);
558       local_size += subset_size;
559     }
560     ierr = PetscFree2(dummy_idx,work);CHKERRQ(ierr);
561     /* free */
562     ierr = ISDestroy(&is_I);CHKERRQ(ierr);
563     ierr = MatDestroy(&AE_II);CHKERRQ(ierr);
564     ierr = PetscFree(all_local_idx_N);CHKERRQ(ierr);
565   } else {
566     Mat         A,F;
567     IS          is_A_all;
568     PetscScalar *work;
569     PetscInt    *idxs_schur,n_I,n_I_all,*dummy_idx;
570 
571     /* get working mat */
572     ierr = ISGetLocalSize(is_I_layer,&n_I);CHKERRQ(ierr);
573     ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size+n_I,all_local_idx_N,PETSC_COPY_VALUES,&is_A_all);CHKERRQ(ierr);
574     ierr = MatGetSubMatrixUnsorted(sub_schurs->A,is_A_all,is_A_all,&A);CHKERRQ(ierr);
575     ierr = ISDestroy(&is_A_all);CHKERRQ(ierr);
576     ierr = MatSetOptionsPrefix(A,"sub_schurs_");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 = ISDestroy(&is_I_layer);CHKERRQ(ierr);
742   ierr = PetscFree(nnz);CHKERRQ(ierr);
743   ierr = MatDestroy(&S_all);CHKERRQ(ierr);
744   ierr = MatDestroy(&S_all_inv);CHKERRQ(ierr);
745   ierr = MatDestroy(&A_BB);CHKERRQ(ierr);
746   ierr = MatDestroy(&A_IB);CHKERRQ(ierr);
747   ierr = MatDestroy(&A_BI);CHKERRQ(ierr);
748   ierr = MatAssemblyBegin(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
749   ierr = MatAssemblyEnd(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
750   if (compute_Stilda) {
751     ierr = MatAssemblyBegin(S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
752     ierr = MatAssemblyEnd(S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
753     ierr = MatAssemblyBegin(S_Ej_inv_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
754     ierr = MatAssemblyEnd(S_Ej_inv_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
755   }
756 
757   /* Global matrix of all assembled Schur on subsets */
758   ierr = MatISSetLocalMat(work_mat,sub_schurs->S_Ej_all);CHKERRQ(ierr);
759   ierr = MatISSetMPIXAIJPreallocation_Private(work_mat,global_schur_subsets,PETSC_TRUE);CHKERRQ(ierr);
760   ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
761 
762   /* Get local part of (\sum_j S_Ej) */
763   ierr = ISCreateGeneral(comm_n,local_size,all_local_idx_G,PETSC_OWN_POINTER,&temp_is);CHKERRQ(ierr);
764   ierr = MatDestroy(&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
765   ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
766 
767   /* Compute explicitly (\sum_j S_Ej)^-1 (faster scaling during PCApply, needs extra work when doing setup) */
768   if (faster_deluxe) {
769     Mat         tmpmat;
770     PetscScalar *array;
771     PetscInt    cum;
772 
773     ierr = MatSeqAIJGetArray(sub_schurs->sum_S_Ej_all,&array);CHKERRQ(ierr);
774     cum = 0;
775     for (i=0;i<sub_schurs->n_subs;i++) {
776       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
777       ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr);
778       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
779       if (!sub_schurs->is_hermitian) {
780         PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,array+cum,&B_N,pivots,&B_ierr));
781         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
782         PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,array+cum,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
783         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
784       } else {
785         PetscInt j,k;
786 
787         PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,array+cum,&B_N,&B_ierr));
788         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
789         PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,array+cum,&B_N,&B_ierr));
790         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
791         for (j=0;j<B_N;j++) {
792           for (k=j+1;k<B_N;k++) {
793             array[k*B_N+j+cum] = array[j*B_N+k+cum];
794           }
795         }
796       }
797       ierr = PetscFPTrapPop();CHKERRQ(ierr);
798       cum += subset_size*subset_size;
799     }
800     ierr = MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_all,&array);CHKERRQ(ierr);
801     ierr = MatMatMult(sub_schurs->S_Ej_all,sub_schurs->sum_S_Ej_all,MAT_INITIAL_MATRIX,1.0,&tmpmat);CHKERRQ(ierr);
802     ierr = MatDestroy(&sub_schurs->S_Ej_all);CHKERRQ(ierr);
803     ierr = MatDestroy(&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
804     sub_schurs->S_Ej_all = tmpmat;
805   }
806 
807   /* Get local part of (\sum_j S^-1_Ej) (\sum_j St^-1_Ej) */
808   if (compute_Stilda) {
809     ierr = MatISSetLocalMat(work_mat,S_Ej_tilda_all);CHKERRQ(ierr);
810     ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
811     ierr = MatDestroy(&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
812     ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
813     ierr = MatISSetLocalMat(work_mat,S_Ej_inv_all);CHKERRQ(ierr);
814     ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
815     ierr = MatDestroy(&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
816     ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
817   }
818 
819   /* free workspace */
820   ierr = PetscFree2(Bwork,pivots);CHKERRQ(ierr);
821   ierr = MatDestroy(&global_schur_subsets);CHKERRQ(ierr);
822   ierr = MatDestroy(&S_Ej_tilda_all);CHKERRQ(ierr);
823   ierr = MatDestroy(&S_Ej_inv_all);CHKERRQ(ierr);
824   ierr = MatDestroy(&work_mat);CHKERRQ(ierr);
825   ierr = ISDestroy(&temp_is);CHKERRQ(ierr);
826   ierr = PetscCommDestroy(&comm_n);CHKERRQ(ierr);
827   PetscFunctionReturn(0);
828 }
829 
830 #undef __FUNCT__
831 #define __FUNCT__ "PCBDDCSubSchursInit"
832 PetscErrorCode PCBDDCSubSchursInit(PCBDDCSubSchurs sub_schurs, IS is_I, IS is_B, PCBDDCGraph graph, ISLocalToGlobalMapping BtoNmap)
833 {
834   IS              *faces,*edges,*all_cc,vertices;
835   PetscInt        i,n_faces,n_edges,n_all_cc;
836   PetscBool       is_sorted;
837   PetscErrorCode  ierr;
838 
839   PetscFunctionBegin;
840   ierr = ISSorted(is_I,&is_sorted);CHKERRQ(ierr);
841   if (!is_sorted) {
842     SETERRQ(PetscObjectComm((PetscObject)is_I),PETSC_ERR_PLIB,"IS for I dofs should be shorted");
843   }
844   ierr = ISSorted(is_B,&is_sorted);CHKERRQ(ierr);
845   if (!is_sorted) {
846     SETERRQ(PetscObjectComm((PetscObject)is_B),PETSC_ERR_PLIB,"IS for B dofs should be shorted");
847   }
848 
849   /* reset any previous data */
850   ierr = PCBDDCSubSchursReset(sub_schurs);CHKERRQ(ierr);
851 
852   /* get index sets for faces and edges (already sorted by global ordering) */
853   ierr = PCBDDCGraphGetCandidatesIS(graph,&n_faces,&faces,&n_edges,&edges,&vertices);CHKERRQ(ierr);
854   n_all_cc = n_faces+n_edges;
855   ierr = PetscBTCreate(n_all_cc,&sub_schurs->is_edge);CHKERRQ(ierr);
856   ierr = PetscBTCreate(n_all_cc,&sub_schurs->computed_Stilda_subs);CHKERRQ(ierr);
857   ierr = PetscMalloc1(n_all_cc,&all_cc);CHKERRQ(ierr);
858   for (i=0;i<n_faces;i++) {
859     all_cc[i] = faces[i];
860   }
861   for (i=0;i<n_edges;i++) {
862     all_cc[n_faces+i] = edges[i];
863     ierr = PetscBTSet(sub_schurs->is_edge,n_faces+i);CHKERRQ(ierr);
864   }
865   ierr = PetscFree(faces);CHKERRQ(ierr);
866   ierr = PetscFree(edges);CHKERRQ(ierr);
867 
868   /* Determine if MUMPS cane be used */
869   sub_schurs->use_mumps = PETSC_FALSE;
870 #if defined(PETSC_HAVE_MUMPS)
871   sub_schurs->use_mumps = PETSC_TRUE;
872 #endif
873 
874   ierr = PetscObjectReference((PetscObject)is_I);CHKERRQ(ierr);
875   sub_schurs->is_I = is_I;
876   ierr = PetscObjectReference((PetscObject)is_B);CHKERRQ(ierr);
877   sub_schurs->is_B = is_B;
878   ierr = PetscObjectReference((PetscObject)graph->l2gmap);CHKERRQ(ierr);
879   sub_schurs->l2gmap = graph->l2gmap;
880   ierr = PetscObjectReference((PetscObject)BtoNmap);CHKERRQ(ierr);
881   sub_schurs->BtoNmap = BtoNmap;
882   sub_schurs->n_subs = n_all_cc;
883   sub_schurs->is_subs = all_cc;
884   if (!sub_schurs->use_mumps) { /* sort by local ordering mumps is not present */
885     for (i=0;i<sub_schurs->n_subs;i++) {
886       ierr = ISSort(sub_schurs->is_subs[i]);CHKERRQ(ierr);
887     }
888   }
889   sub_schurs->is_Ej_com = vertices;
890 
891   /* allocate space for schur complements */
892   sub_schurs->S_Ej_all = NULL;
893   sub_schurs->sum_S_Ej_all = NULL;
894   sub_schurs->sum_S_Ej_inv_all = NULL;
895   sub_schurs->sum_S_Ej_tilda_all = NULL;
896   sub_schurs->is_Ej_all = NULL;
897   PetscFunctionReturn(0);
898 }
899 
900 #undef __FUNCT__
901 #define __FUNCT__ "PCBDDCSubSchursCreate"
902 PetscErrorCode PCBDDCSubSchursCreate(PCBDDCSubSchurs *sub_schurs)
903 {
904   PCBDDCSubSchurs schurs_ctx;
905   PetscErrorCode  ierr;
906 
907   PetscFunctionBegin;
908   ierr = PetscNew(&schurs_ctx);CHKERRQ(ierr);
909   schurs_ctx->n_subs = 0;
910   *sub_schurs = schurs_ctx;
911   PetscFunctionReturn(0);
912 }
913 
914 #undef __FUNCT__
915 #define __FUNCT__ "PCBDDCSubSchursDestroy"
916 PetscErrorCode PCBDDCSubSchursDestroy(PCBDDCSubSchurs *sub_schurs)
917 {
918   PetscErrorCode ierr;
919 
920   PetscFunctionBegin;
921   ierr = PCBDDCSubSchursReset(*sub_schurs);CHKERRQ(ierr);
922   ierr = PetscFree(*sub_schurs);CHKERRQ(ierr);
923   PetscFunctionReturn(0);
924 }
925 
926 #undef __FUNCT__
927 #define __FUNCT__ "PCBDDCSubSchursReset"
928 PetscErrorCode PCBDDCSubSchursReset(PCBDDCSubSchurs sub_schurs)
929 {
930   PetscInt       i;
931   PetscErrorCode ierr;
932 
933   PetscFunctionBegin;
934   ierr = MatDestroy(&sub_schurs->A);CHKERRQ(ierr);
935   ierr = MatDestroy(&sub_schurs->S);CHKERRQ(ierr);
936   ierr = ISDestroy(&sub_schurs->is_I);CHKERRQ(ierr);
937   ierr = ISDestroy(&sub_schurs->is_B);CHKERRQ(ierr);
938   ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->l2gmap);CHKERRQ(ierr);
939   ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->BtoNmap);CHKERRQ(ierr);
940   ierr = MatDestroy(&sub_schurs->S_Ej_all);CHKERRQ(ierr);
941   ierr = MatDestroy(&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
942   ierr = MatDestroy(&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
943   ierr = MatDestroy(&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
944   ierr = ISDestroy(&sub_schurs->is_Ej_all);CHKERRQ(ierr);
945   ierr = ISDestroy(&sub_schurs->is_Ej_com);CHKERRQ(ierr);
946   ierr = PetscBTDestroy(&sub_schurs->is_edge);CHKERRQ(ierr);
947   ierr = PetscBTDestroy(&sub_schurs->computed_Stilda_subs);CHKERRQ(ierr);
948   for (i=0;i<sub_schurs->n_subs;i++) {
949     ierr = ISDestroy(&sub_schurs->is_subs[i]);CHKERRQ(ierr);
950   }
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