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