xref: /petsc/src/ksp/pc/impls/bddc/bddcschurs.c (revision 6dba178d53d088f8af27c39b86adaa0496f0d2e2)
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__ "PCBDDCMumpsReuseSolve_Private"
12 static PetscErrorCode PCBDDCMumpsReuseSolve_Private(PC pc, Vec rhs, Vec sol, PetscBool transpose, PetscBool full)
13 {
14   PCBDDCReuseMumps ctx;
15   PetscScalar      *array,*array_mumps;
16   PetscInt         n;
17   PetscBool        copy = PETSC_FALSE;
18   PetscErrorCode   ierr;
19 
20   PetscFunctionBegin;
21   ierr = PCShellGetContext(pc,(void **)&ctx);CHKERRQ(ierr);
22   if (full) {
23 #if defined(PETSC_HAVE_MUMPS)
24     ierr = MatMumpsSetIcntl(ctx->F,26,-1);CHKERRQ(ierr);
25 #endif
26     copy = ctx->has_vertices;
27   } else {
28 #if defined(PETSC_HAVE_MUMPS)
29     ierr = MatMumpsSetIcntl(ctx->F,26,0);CHKERRQ(ierr);
30 #endif
31     copy = PETSC_TRUE;
32   }
33   ierr = VecGetLocalSize(rhs,&n);CHKERRQ(ierr);
34   /* copy rhs into factored matrix workspace */
35   if (copy) {
36     ierr = VecGetArrayRead(rhs,(const PetscScalar**)&array);CHKERRQ(ierr);
37     ierr = VecGetArray(ctx->rhs,&array_mumps);CHKERRQ(ierr);
38     ierr = PetscMemcpy(array_mumps,array,n*sizeof(PetscScalar));CHKERRQ(ierr);
39     ierr = VecRestoreArray(ctx->rhs,&array_mumps);CHKERRQ(ierr);
40     ierr = VecRestoreArrayRead(rhs,(const PetscScalar**)&array);CHKERRQ(ierr);
41 
42     if (transpose) {
43       ierr = MatSolveTranspose(ctx->F,ctx->rhs,ctx->sol);CHKERRQ(ierr);
44     } else {
45       ierr = MatSolve(ctx->F,ctx->rhs,ctx->sol);CHKERRQ(ierr);
46     }
47 
48     /* get back data to caller worskpace */
49     ierr = VecGetArrayRead(ctx->sol,(const PetscScalar**)&array_mumps);CHKERRQ(ierr);
50     ierr = VecGetArray(sol,&array);CHKERRQ(ierr);
51     ierr = PetscMemcpy(array,array_mumps,n*sizeof(PetscScalar));CHKERRQ(ierr);
52     ierr = VecRestoreArray(sol,&array);CHKERRQ(ierr);
53     ierr = VecRestoreArrayRead(ctx->sol,(const PetscScalar**)&array_mumps);CHKERRQ(ierr);
54   } else {
55     if (transpose) {
56       ierr = MatSolveTranspose(ctx->F,rhs,sol);CHKERRQ(ierr);
57     } else {
58       ierr = MatSolve(ctx->F,rhs,sol);CHKERRQ(ierr);
59     }
60   }
61   PetscFunctionReturn(0);
62 }
63 
64 #undef __FUNCT__
65 #define __FUNCT__ "PCBDDCMumpsCorrectionSolve"
66 static PetscErrorCode PCBDDCMumpsCorrectionSolve(PC pc, Vec rhs, Vec sol)
67 {
68   PetscErrorCode   ierr;
69 
70   PetscFunctionBegin;
71   ierr = PCBDDCMumpsReuseSolve_Private(pc,rhs,sol,PETSC_FALSE,PETSC_TRUE);CHKERRQ(ierr);
72   PetscFunctionReturn(0);
73 }
74 
75 #undef __FUNCT__
76 #define __FUNCT__ "PCBDDCMumpsCorrectionSolveTranspose"
77 static PetscErrorCode PCBDDCMumpsCorrectionSolveTranspose(PC pc, Vec rhs, Vec sol)
78 {
79   PetscErrorCode   ierr;
80 
81   PetscFunctionBegin;
82   ierr = PCBDDCMumpsReuseSolve_Private(pc,rhs,sol,PETSC_TRUE,PETSC_TRUE);CHKERRQ(ierr);
83   PetscFunctionReturn(0);
84 }
85 
86 #undef __FUNCT__
87 #define __FUNCT__ "PCBDDCMumpsInteriorSolve"
88 static PetscErrorCode PCBDDCMumpsInteriorSolve(PC pc, Vec rhs, Vec sol)
89 {
90   PetscErrorCode   ierr;
91 
92   PetscFunctionBegin;
93   ierr = PCBDDCMumpsReuseSolve_Private(pc,rhs,sol,PETSC_FALSE,PETSC_FALSE);CHKERRQ(ierr);
94   PetscFunctionReturn(0);
95 }
96 
97 #undef __FUNCT__
98 #define __FUNCT__ "PCBDDCMumpsInteriorSolveTranspose"
99 static PetscErrorCode PCBDDCMumpsInteriorSolveTranspose(PC pc, Vec rhs, Vec sol)
100 {
101   PetscErrorCode   ierr;
102 
103   PetscFunctionBegin;
104   ierr = PCBDDCMumpsReuseSolve_Private(pc,rhs,sol,PETSC_TRUE,PETSC_FALSE);CHKERRQ(ierr);
105   PetscFunctionReturn(0);
106 }
107 
108 #undef __FUNCT__
109 #define __FUNCT__ "PCBDDCReuseMumpsReset"
110 static PetscErrorCode PCBDDCReuseMumpsReset(PCBDDCReuseMumps reuse)
111 {
112   PetscErrorCode ierr;
113 
114   PetscFunctionBegin;
115   ierr = MatDestroy(&reuse->F);CHKERRQ(ierr);
116   ierr = VecDestroy(&reuse->sol);CHKERRQ(ierr);
117   ierr = VecDestroy(&reuse->rhs);CHKERRQ(ierr);
118   ierr = PCDestroy(&reuse->interior_solver);CHKERRQ(ierr);
119   ierr = PCDestroy(&reuse->correction_solver);CHKERRQ(ierr);
120   ierr = ISDestroy(&reuse->is_R);CHKERRQ(ierr);
121   ierr = ISDestroy(&reuse->is_B);CHKERRQ(ierr);
122   ierr = VecScatterDestroy(&reuse->correction_scatter_B);CHKERRQ(ierr);
123   ierr = VecDestroy(&reuse->sol_B);CHKERRQ(ierr);
124   ierr = VecDestroy(&reuse->rhs_B);CHKERRQ(ierr);
125   PetscFunctionReturn(0);
126 }
127 
128 #undef __FUNCT__
129 #define __FUNCT__ "PCBDDCComputeExplicitSchur"
130 static PetscErrorCode PCBDDCComputeExplicitSchur(Mat M, PetscBool issym, MatReuse reuse, Mat *S)
131 {
132   Mat            B, C, D, Bd, Cd, AinvBd;
133   KSP            ksp;
134   PC             pc;
135   PetscBool      isLU, isILU, isCHOL, Bdense, Cdense;
136   PetscReal      fill = 2.0;
137   PetscInt       n_I;
138   PetscMPIInt    size;
139   PetscErrorCode ierr;
140 
141   PetscFunctionBegin;
142   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)M),&size);CHKERRQ(ierr);
143   if (size != 1) {
144     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not for parallel matrices");
145   }
146   if (reuse == MAT_REUSE_MATRIX) {
147     PetscBool Sdense;
148 
149     ierr = PetscObjectTypeCompare((PetscObject)*S, MATSEQDENSE, &Sdense);CHKERRQ(ierr);
150     if (!Sdense) {
151       SETERRQ(PetscObjectComm((PetscObject)M),PETSC_ERR_SUP,"S should be dense");
152     }
153   }
154   ierr = MatSchurComplementGetSubMatrices(M, NULL, NULL, &B, &C, &D);CHKERRQ(ierr);
155   ierr = MatSchurComplementGetKSP(M, &ksp);CHKERRQ(ierr);
156   ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr);
157   ierr = PetscObjectTypeCompare((PetscObject) pc, PCLU, &isLU);CHKERRQ(ierr);
158   ierr = PetscObjectTypeCompare((PetscObject) pc, PCILU, &isILU);CHKERRQ(ierr);
159   ierr = PetscObjectTypeCompare((PetscObject) pc, PCCHOLESKY, &isCHOL);CHKERRQ(ierr);
160   ierr = PetscObjectTypeCompare((PetscObject) B, MATSEQDENSE, &Bdense);CHKERRQ(ierr);
161   ierr = PetscObjectTypeCompare((PetscObject) C, MATSEQDENSE, &Cdense);CHKERRQ(ierr);
162   ierr = MatGetSize(B,&n_I,NULL);CHKERRQ(ierr);
163   if (n_I) {
164     if (!Bdense) {
165       ierr = MatConvert(B, MATSEQDENSE, MAT_INITIAL_MATRIX, &Bd);CHKERRQ(ierr);
166     } else {
167       Bd = B;
168     }
169 
170     if (isLU || isILU || isCHOL) {
171       Mat fact;
172       ierr = KSPSetUp(ksp);CHKERRQ(ierr);
173       ierr = PCFactorGetMatrix(pc, &fact);CHKERRQ(ierr);
174       ierr = MatDuplicate(Bd, MAT_DO_NOT_COPY_VALUES, &AinvBd);CHKERRQ(ierr);
175       ierr = MatMatSolve(fact, Bd, AinvBd);CHKERRQ(ierr);
176     } else {
177       PetscBool ex = PETSC_TRUE;
178 
179       if (ex) {
180         Mat Ainvd;
181 
182         ierr = PCComputeExplicitOperator(pc, &Ainvd);CHKERRQ(ierr);
183         ierr = MatMatMult(Ainvd, Bd, MAT_INITIAL_MATRIX, fill, &AinvBd);CHKERRQ(ierr);
184         ierr = MatDestroy(&Ainvd);CHKERRQ(ierr);
185       } else {
186         Vec         sol,rhs;
187         PetscScalar *arrayrhs,*arraysol;
188         PetscInt    i,nrhs,n;
189 
190         ierr = MatDuplicate(Bd, MAT_DO_NOT_COPY_VALUES, &AinvBd);CHKERRQ(ierr);
191         ierr = MatGetSize(Bd,&n,&nrhs);CHKERRQ(ierr);
192         ierr = MatDenseGetArray(Bd,&arrayrhs);CHKERRQ(ierr);
193         ierr = MatDenseGetArray(AinvBd,&arraysol);CHKERRQ(ierr);
194         ierr = KSPGetSolution(ksp,&sol);CHKERRQ(ierr);
195         ierr = KSPGetRhs(ksp,&rhs);CHKERRQ(ierr);
196         for (i=0;i<nrhs;i++) {
197           ierr = VecPlaceArray(rhs,arrayrhs+i*n);CHKERRQ(ierr);
198           ierr = VecPlaceArray(sol,arraysol+i*n);CHKERRQ(ierr);
199           ierr = KSPSolve(ksp,rhs,sol);CHKERRQ(ierr);
200           ierr = VecResetArray(rhs);CHKERRQ(ierr);
201           ierr = VecResetArray(sol);CHKERRQ(ierr);
202         }
203         ierr = MatDenseRestoreArray(Bd,&arrayrhs);CHKERRQ(ierr);
204         ierr = MatDenseRestoreArray(AinvBd,&arrayrhs);CHKERRQ(ierr);
205       }
206     }
207     if (!Bdense & !issym) {
208       ierr = MatDestroy(&Bd);CHKERRQ(ierr);
209     }
210 
211     if (!issym) {
212       if (!Cdense) {
213         ierr = MatConvert(C, MATSEQDENSE, MAT_INITIAL_MATRIX, &Cd);CHKERRQ(ierr);
214       } else {
215         Cd = C;
216       }
217       ierr = MatMatMult(Cd, AinvBd, reuse, fill, S);CHKERRQ(ierr);
218       if (!Cdense) {
219         ierr = MatDestroy(&Cd);CHKERRQ(ierr);
220       }
221     } else {
222       ierr = MatTransposeMatMult(Bd, AinvBd, reuse, fill, S);CHKERRQ(ierr);
223       if (!Bdense) {
224         ierr = MatDestroy(&Bd);CHKERRQ(ierr);
225       }
226     }
227     ierr = MatDestroy(&AinvBd);CHKERRQ(ierr);
228   }
229 
230   if (D) {
231     Mat       Dd;
232     PetscBool Ddense;
233 
234     ierr = PetscObjectTypeCompare((PetscObject)D,MATSEQDENSE,&Ddense);CHKERRQ(ierr);
235     if (!Ddense) {
236       ierr = MatConvert(D, MATSEQDENSE, MAT_INITIAL_MATRIX, &Dd);CHKERRQ(ierr);
237     } else {
238       Dd = D;
239     }
240     if (n_I) {
241       ierr = MatAYPX(*S,-1.0,Dd,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
242     } else {
243       if (reuse == MAT_INITIAL_MATRIX) {
244         ierr = MatDuplicate(Dd,MAT_COPY_VALUES,S);CHKERRQ(ierr);
245       } else {
246         ierr = MatCopy(Dd,*S,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
247       }
248     }
249     if (!Ddense) {
250       ierr = MatDestroy(&Dd);CHKERRQ(ierr);
251     }
252   } else {
253     ierr = MatScale(*S,-1.0);CHKERRQ(ierr);
254   }
255   PetscFunctionReturn(0);
256 }
257 
258 #undef __FUNCT__
259 #define __FUNCT__ "PCBDDCSubSchursSetUp"
260 PetscErrorCode PCBDDCSubSchursSetUp(PCBDDCSubSchurs sub_schurs, Mat Ain, Mat Sin, PetscBool exact_schur, PetscInt xadj[], PetscInt adjncy[], PetscInt nlayers, PetscBool faster_deluxe, PetscBool compute_Stilda, PetscBool reuse_solvers, PetscBool benign_trick)
261 {
262   Mat                    F,A_II,A_IB,A_BI,A_BB,AE_II;
263   Mat                    S_all;
264   Mat                    global_schur_subsets,work_mat,*submats;
265   ISLocalToGlobalMapping l2gmap_subsets;
266   IS                     is_I,is_I_layer;
267   IS                     all_subsets,all_subsets_mult,all_subsets_n;
268   PetscInt               *nnz,*all_local_idx_N;
269   PetscInt               *auxnum1,*auxnum2;
270   PetscInt               i,subset_size,max_subset_size;
271   PetscInt               n_B,extra,local_size,global_size;
272   PetscBLASInt           B_N,B_ierr,B_lwork,*pivots;
273   PetscScalar            *Bwork;
274   PetscSubcomm           subcomm;
275   PetscMPIInt            color,rank;
276   MPI_Comm               comm_n;
277   PetscBool              use_getr = PETSC_FALSE;
278   PetscErrorCode         ierr;
279 
280   PetscFunctionBegin;
281   /* update info in sub_schurs */
282   ierr = MatDestroy(&sub_schurs->A);CHKERRQ(ierr);
283   ierr = MatDestroy(&sub_schurs->S);CHKERRQ(ierr);
284   sub_schurs->is_hermitian = PETSC_FALSE;
285   sub_schurs->is_posdef = PETSC_FALSE;
286   if (Ain) {
287     PetscBool isseqaij;
288     /* determine if we are dealing with hermitian positive definite problems */
289 #if !defined(PETSC_USE_COMPLEX)
290     if (Ain->symmetric_set) {
291       sub_schurs->is_hermitian = Ain->symmetric;
292     }
293 #else
294     if (Ain->hermitian_set) {
295       sub_schurs->is_hermitian = Ain->hermitian;
296     }
297 #endif
298     if (Ain->spd_set) {
299       sub_schurs->is_posdef = Ain->spd;
300     }
301 
302     /* check */
303     ierr = PetscObjectTypeCompare((PetscObject)Ain,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
304     if (compute_Stilda && (!Ain->hermitian_set || !Ain->spd_set)) { /* these are lazy checks, maybe I should throw an error if they are not set */
305       PetscInt lsize;
306 
307       ierr = MatGetSize(Ain,&lsize,NULL);CHKERRQ(ierr);
308       if (lsize) {
309         PetscScalar val;
310         PetscReal   norm;
311         Vec         vec1,vec2,vec3;
312 
313         ierr = MatCreateVecs(Ain,&vec1,&vec2);CHKERRQ(ierr);
314         ierr = VecDuplicate(vec1,&vec3);CHKERRQ(ierr);
315         ierr = VecSetRandom(vec1,NULL);CHKERRQ(ierr);
316         ierr = MatMult(Ain,vec1,vec2);CHKERRQ(ierr);
317 #if !defined(PETSC_USE_COMPLEX)
318         ierr = MatMultTranspose(Ain,vec1,vec3);CHKERRQ(ierr);
319 #else
320         ierr = MatMultHermitianTranspose(Ain,vec1,vec3);CHKERRQ(ierr);
321 #endif
322         ierr = VecAXPY(vec3,-1.0,vec2);CHKERRQ(ierr);
323         ierr = VecNorm(vec3,NORM_INFINITY,&norm);CHKERRQ(ierr);
324         if (!Ain->hermitian_set) {
325           if (norm > PetscSqrtReal(PETSC_SMALL)) {
326             sub_schurs->is_hermitian = PETSC_FALSE;
327           } else {
328             sub_schurs->is_hermitian = PETSC_TRUE;
329           }
330         }
331         if (!Ain->spd_set && !benign_trick) {
332           ierr = VecDot(vec1,vec2,&val);CHKERRQ(ierr);
333           if (PetscRealPart(val) > 0. && PetscAbsReal(PetscImaginaryPart(val)) < PETSC_SMALL) sub_schurs->is_posdef = PETSC_TRUE;
334         }
335         ierr = VecDestroy(&vec1);CHKERRQ(ierr);
336         ierr = VecDestroy(&vec2);CHKERRQ(ierr);
337         ierr = VecDestroy(&vec3);CHKERRQ(ierr);
338       } else {
339         sub_schurs->is_hermitian = PETSC_TRUE;
340         sub_schurs->is_posdef = PETSC_TRUE;
341       }
342     }
343     if (isseqaij) {
344       ierr = PetscObjectReference((PetscObject)Ain);CHKERRQ(ierr);
345       sub_schurs->A = Ain;
346     } else {
347       ierr = MatConvert(Ain,MATSEQAIJ,MAT_INITIAL_MATRIX,&sub_schurs->A);CHKERRQ(ierr);
348     }
349   }
350 
351   ierr = PetscObjectReference((PetscObject)Sin);CHKERRQ(ierr);
352   sub_schurs->S = Sin;
353   if (sub_schurs->use_mumps) {
354     sub_schurs->use_mumps = (PetscBool)(!!sub_schurs->A);
355   }
356 
357   /* preliminary checks */
358   if (!sub_schurs->use_mumps && compute_Stilda) {
359     SETERRQ(PetscObjectComm((PetscObject)sub_schurs->l2gmap),PETSC_ERR_SUP,"Adaptive selection of constraints requires MUMPS");
360   }
361 
362   /* restrict work on active processes */
363   color = 0;
364   if (!sub_schurs->n_subs) color = 1; /* this can happen if we are in a multilevel case or if the subdomain is disconnected */
365   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&rank);CHKERRQ(ierr);
366   ierr = PetscSubcommCreate(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&subcomm);CHKERRQ(ierr);
367   ierr = PetscSubcommSetNumber(subcomm,2);CHKERRQ(ierr);
368   ierr = PetscSubcommSetTypeGeneral(subcomm,color,rank);CHKERRQ(ierr);
369   ierr = PetscCommDuplicate(PetscSubcommChild(subcomm),&comm_n,NULL);CHKERRQ(ierr);
370   ierr = PetscSubcommDestroy(&subcomm);CHKERRQ(ierr);
371   if (!sub_schurs->n_subs) {
372     ierr = PetscCommDestroy(&comm_n);CHKERRQ(ierr);
373     PetscFunctionReturn(0);
374   }
375   /* ierr = PetscCommDuplicate(PetscObjectComm((PetscObject)sub_schurs->l2gmap),&comm_n,NULL);CHKERRQ(ierr); */
376 
377   /* get Schur complement matrices */
378   if (!sub_schurs->use_mumps) {
379     Mat       tA_IB,tA_BI,tA_BB;
380     PetscBool isseqsbaij;
381     ierr = MatSchurComplementGetSubMatrices(sub_schurs->S,&A_II,NULL,&tA_IB,&tA_BI,&tA_BB);CHKERRQ(ierr);
382     ierr = PetscObjectTypeCompare((PetscObject)tA_BB,MATSEQSBAIJ,&isseqsbaij);CHKERRQ(ierr);
383     if (isseqsbaij) {
384       ierr = MatConvert(tA_BB,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_BB);CHKERRQ(ierr);
385       ierr = MatConvert(tA_IB,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_IB);CHKERRQ(ierr);
386       ierr = MatConvert(tA_BI,MATSEQAIJ,MAT_INITIAL_MATRIX,&A_BI);CHKERRQ(ierr);
387     } else {
388       ierr = PetscObjectReference((PetscObject)tA_BB);CHKERRQ(ierr);
389       A_BB = tA_BB;
390       ierr = PetscObjectReference((PetscObject)tA_IB);CHKERRQ(ierr);
391       A_IB = tA_IB;
392       ierr = PetscObjectReference((PetscObject)tA_BI);CHKERRQ(ierr);
393       A_BI = tA_BI;
394     }
395   } else {
396     A_II = NULL;
397     A_IB = NULL;
398     A_BI = NULL;
399     A_BB = NULL;
400   }
401   S_all = NULL;
402 
403   /* determine interior problems */
404   ierr = ISGetLocalSize(sub_schurs->is_I,&i);CHKERRQ(ierr);
405   if (nlayers >= 0 && i) { /* Interior problems can be different from the original one */
406     PetscBT                touched;
407     const PetscInt*        idx_B;
408     PetscInt               n_I,n_B,n_local_dofs,n_prev_added,j,layer,*local_numbering;
409 
410     if (xadj == NULL || adjncy == NULL) {
411       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cannot request layering without adjacency");
412     }
413     /* get sizes */
414     ierr = ISGetLocalSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr);
415     ierr = ISGetLocalSize(sub_schurs->is_B,&n_B);CHKERRQ(ierr);
416 
417     ierr = PetscMalloc1(n_I+n_B,&local_numbering);CHKERRQ(ierr);
418     ierr = PetscBTCreate(n_I+n_B,&touched);CHKERRQ(ierr);
419     ierr = PetscBTMemzero(n_I+n_B,touched);CHKERRQ(ierr);
420 
421     /* all boundary dofs must be skipped when adding layers */
422     ierr = ISGetIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr);
423     for (j=0;j<n_B;j++) {
424       ierr = PetscBTSet(touched,idx_B[j]);CHKERRQ(ierr);
425     }
426     ierr = PetscMemcpy(local_numbering,idx_B,n_B*sizeof(PetscInt));CHKERRQ(ierr);
427     ierr = ISRestoreIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr);
428 
429     /* add prescribed number of layers of dofs */
430     n_local_dofs = n_B;
431     n_prev_added = n_B;
432     for (layer=0;layer<nlayers;layer++) {
433       PetscInt n_added;
434       if (n_local_dofs == n_I+n_B) break;
435       if (n_local_dofs > n_I+n_B) {
436         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);
437       }
438       ierr = PCBDDCAdjGetNextLayer_Private(local_numbering+n_local_dofs,n_prev_added,touched,xadj,adjncy,&n_added);CHKERRQ(ierr);
439       n_prev_added = n_added;
440       n_local_dofs += n_added;
441       if (!n_added) break;
442     }
443     ierr = PetscBTDestroy(&touched);CHKERRQ(ierr);
444 
445     /* IS for I layer dofs in original numbering */
446     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);
447     ierr = PetscFree(local_numbering);CHKERRQ(ierr);
448     ierr = ISSort(is_I_layer);CHKERRQ(ierr);
449     /* IS for I layer dofs in I numbering */
450     if (!sub_schurs->use_mumps) {
451       ISLocalToGlobalMapping ItoNmap;
452       ierr = ISLocalToGlobalMappingCreateIS(sub_schurs->is_I,&ItoNmap);CHKERRQ(ierr);
453       ierr = ISGlobalToLocalMappingApplyIS(ItoNmap,IS_GTOLM_DROP,is_I_layer,&is_I);CHKERRQ(ierr);
454       ierr = ISLocalToGlobalMappingDestroy(&ItoNmap);CHKERRQ(ierr);
455 
456       /* II block */
457       ierr = MatGetSubMatrix(A_II,is_I,is_I,MAT_INITIAL_MATRIX,&AE_II);CHKERRQ(ierr);
458     }
459   } else {
460     PetscInt n_I;
461 
462     /* IS for I dofs in original numbering */
463     ierr = PetscObjectReference((PetscObject)sub_schurs->is_I);CHKERRQ(ierr);
464     is_I_layer = sub_schurs->is_I;
465 
466     /* IS for I dofs in I numbering (strided 1) */
467     if (!sub_schurs->use_mumps) {
468       ierr = ISGetSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr);
469       ierr = ISCreateStride(PetscObjectComm((PetscObject)sub_schurs->is_I),n_I,0,1,&is_I);CHKERRQ(ierr);
470 
471       /* II block is the same */
472       ierr = PetscObjectReference((PetscObject)A_II);CHKERRQ(ierr);
473       AE_II = A_II;
474     }
475   }
476 
477   /* Get info on subset sizes and sum of all subsets sizes */
478   max_subset_size = 0;
479   local_size = 0;
480   for (i=0;i<sub_schurs->n_subs;i++) {
481     ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
482     max_subset_size = PetscMax(subset_size,max_subset_size);
483     local_size += subset_size;
484   }
485 
486   /* Work arrays for local indices */
487   extra = 0;
488   ierr = ISGetLocalSize(sub_schurs->is_B,&n_B);CHKERRQ(ierr);
489   if (sub_schurs->use_mumps && is_I_layer) {
490     ierr = ISGetLocalSize(is_I_layer,&extra);CHKERRQ(ierr);
491   }
492   ierr = PetscMalloc1(n_B+extra,&all_local_idx_N);CHKERRQ(ierr);
493   if (extra) {
494     const PetscInt *idxs;
495     ierr = ISGetIndices(is_I_layer,&idxs);CHKERRQ(ierr);
496     ierr = PetscMemcpy(all_local_idx_N,idxs,extra*sizeof(PetscInt));CHKERRQ(ierr);
497     ierr = ISRestoreIndices(is_I_layer,&idxs);CHKERRQ(ierr);
498   }
499   ierr = PetscMalloc1(local_size,&nnz);CHKERRQ(ierr);
500   ierr = PetscMalloc1(sub_schurs->n_subs,&auxnum1);CHKERRQ(ierr);
501   ierr = PetscMalloc1(sub_schurs->n_subs,&auxnum2);CHKERRQ(ierr);
502 
503   /* Get local indices in local numbering */
504   local_size = 0;
505   for (i=0;i<sub_schurs->n_subs;i++) {
506     PetscInt j;
507     const    PetscInt *idxs;
508 
509     ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
510     ierr = ISGetIndices(sub_schurs->is_subs[i],&idxs);CHKERRQ(ierr);
511     /* start (smallest in global ordering) and multiplicity */
512     auxnum1[i] = idxs[0];
513     auxnum2[i] = subset_size;
514     /* subset indices in local numbering */
515     ierr = PetscMemcpy(all_local_idx_N+local_size+extra,idxs,subset_size*sizeof(PetscInt));CHKERRQ(ierr);
516     ierr = ISRestoreIndices(sub_schurs->is_subs[i],&idxs);CHKERRQ(ierr);
517     for (j=0;j<subset_size;j++) nnz[local_size+j] = subset_size;
518     local_size += subset_size;
519   }
520 
521   /* allocate extra workspace needed only for GETRI */
522   Bwork = NULL;
523   pivots = NULL;
524   if (local_size && !benign_trick && (!sub_schurs->is_hermitian || !sub_schurs->is_posdef)) {
525     PetscScalar lwork;
526 
527     use_getr = PETSC_TRUE;
528     B_lwork = -1;
529     ierr = PetscBLASIntCast(local_size,&B_N);CHKERRQ(ierr);
530     ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
531     PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,Bwork,&B_N,pivots,&lwork,&B_lwork,&B_ierr));
532     ierr = PetscFPTrapPop();CHKERRQ(ierr);
533     if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to GETRI Lapack routine %d",(int)B_ierr);
534     ierr = PetscBLASIntCast((PetscInt)PetscRealPart(lwork),&B_lwork);CHKERRQ(ierr);
535     ierr = PetscMalloc2(B_lwork,&Bwork,B_N,&pivots);CHKERRQ(ierr);
536   }
537 
538   /* prepare parallel matrices for summing up properly schurs on subsets */
539   ierr = ISCreateGeneral(comm_n,sub_schurs->n_subs,auxnum1,PETSC_OWN_POINTER,&all_subsets_n);CHKERRQ(ierr);
540   ierr = ISLocalToGlobalMappingApplyIS(sub_schurs->l2gmap,all_subsets_n,&all_subsets);CHKERRQ(ierr);
541   ierr = ISDestroy(&all_subsets_n);CHKERRQ(ierr);
542   ierr = ISCreateGeneral(comm_n,sub_schurs->n_subs,auxnum2,PETSC_OWN_POINTER,&all_subsets_mult);CHKERRQ(ierr);
543   ierr = PCBDDCSubsetNumbering(all_subsets,all_subsets_mult,&global_size,&all_subsets_n);CHKERRQ(ierr);
544   ierr = ISDestroy(&all_subsets);CHKERRQ(ierr);
545   ierr = ISDestroy(&all_subsets_mult);CHKERRQ(ierr);
546   ierr = ISGetLocalSize(all_subsets_n,&i);CHKERRQ(ierr);
547   if (i != local_size) {
548     SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Invalid size of new subset! %d != %d",i,local_size);
549   }
550   ierr = ISLocalToGlobalMappingCreateIS(all_subsets_n,&l2gmap_subsets);CHKERRQ(ierr);
551   ierr = MatCreateIS(comm_n,1,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size,l2gmap_subsets,NULL,&work_mat);CHKERRQ(ierr);
552   ierr = ISLocalToGlobalMappingDestroy(&l2gmap_subsets);CHKERRQ(ierr);
553   ierr = MatCreate(PetscObjectComm((PetscObject)work_mat),&global_schur_subsets);CHKERRQ(ierr);
554   ierr = MatSetSizes(global_schur_subsets,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size);CHKERRQ(ierr);
555   ierr = MatSetType(global_schur_subsets,MATMPIAIJ);CHKERRQ(ierr);
556 
557   /* subset indices in local boundary numbering */
558   if (!sub_schurs->is_Ej_all) {
559     PetscInt *all_local_idx_B;
560 
561     ierr = PetscMalloc1(local_size,&all_local_idx_B);CHKERRQ(ierr);
562     ierr = ISGlobalToLocalMappingApply(sub_schurs->BtoNmap,IS_GTOLM_DROP,local_size,all_local_idx_N+extra,&subset_size,all_local_idx_B);CHKERRQ(ierr);
563     if (subset_size != local_size) {
564       SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in sub_schurs serial (BtoNmap)! %d != %d\n",subset_size,local_size);
565     }
566     ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size,all_local_idx_B,PETSC_OWN_POINTER,&sub_schurs->is_Ej_all);CHKERRQ(ierr);
567   }
568 
569   /* Local matrix of all local Schur on subsets (transposed) */
570   if (!sub_schurs->S_Ej_all) {
571     ierr = MatCreate(PETSC_COMM_SELF,&sub_schurs->S_Ej_all);CHKERRQ(ierr);
572     ierr = MatSetSizes(sub_schurs->S_Ej_all,PETSC_DECIDE,PETSC_DECIDE,local_size,local_size);CHKERRQ(ierr);
573     ierr = MatSetType(sub_schurs->S_Ej_all,MATAIJ);CHKERRQ(ierr);
574     ierr = MatSeqAIJSetPreallocation(sub_schurs->S_Ej_all,0,nnz);CHKERRQ(ierr);
575   }
576 
577   /* Compute Schur complements explicitly */
578   F = NULL;
579   if (!sub_schurs->use_mumps) { /* this code branch is used when MUMPS is not present; it is not very efficient, unless the economic version of the scaling is required */
580     Mat         S_Ej_expl;
581     PetscScalar *work;
582     PetscInt    j,*dummy_idx;
583     PetscBool   Sdense;
584 
585     ierr = PetscMalloc2(max_subset_size,&dummy_idx,max_subset_size*max_subset_size,&work);CHKERRQ(ierr);
586     local_size = 0;
587     for (i=0;i<sub_schurs->n_subs;i++) {
588       IS  is_subset_B;
589       Mat AE_EE,AE_IE,AE_EI,S_Ej;
590 
591       /* subsets in original and boundary numbering */
592       ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,sub_schurs->is_subs[i],&is_subset_B);CHKERRQ(ierr);
593       /* EE block */
594       ierr = MatGetSubMatrix(A_BB,is_subset_B,is_subset_B,MAT_INITIAL_MATRIX,&AE_EE);CHKERRQ(ierr);
595       /* IE block */
596       ierr = MatGetSubMatrix(A_IB,is_I,is_subset_B,MAT_INITIAL_MATRIX,&AE_IE);CHKERRQ(ierr);
597       /* EI block */
598       if (sub_schurs->is_hermitian) {
599         ierr = MatCreateTranspose(AE_IE,&AE_EI);CHKERRQ(ierr);
600       } else {
601         ierr = MatGetSubMatrix(A_BI,is_subset_B,is_I,MAT_INITIAL_MATRIX,&AE_EI);CHKERRQ(ierr);
602       }
603       ierr = ISDestroy(&is_subset_B);CHKERRQ(ierr);
604       ierr = MatCreateSchurComplement(AE_II,AE_II,AE_IE,AE_EI,AE_EE,&S_Ej);CHKERRQ(ierr);
605       ierr = MatDestroy(&AE_EE);CHKERRQ(ierr);
606       ierr = MatDestroy(&AE_IE);CHKERRQ(ierr);
607       ierr = MatDestroy(&AE_EI);CHKERRQ(ierr);
608       if (AE_II == A_II) { /* we can reuse the same ksp */
609         KSP ksp;
610         ierr = MatSchurComplementGetKSP(sub_schurs->S,&ksp);CHKERRQ(ierr);
611         ierr = MatSchurComplementSetKSP(S_Ej,ksp);CHKERRQ(ierr);
612       } else { /* build new ksp object which inherits ksp and pc types from the original one */
613         KSP       origksp,schurksp;
614         PC        origpc,schurpc;
615         KSPType   ksp_type;
616         PetscInt  n_internal;
617         PetscBool ispcnone;
618 
619         ierr = MatSchurComplementGetKSP(sub_schurs->S,&origksp);CHKERRQ(ierr);
620         ierr = MatSchurComplementGetKSP(S_Ej,&schurksp);CHKERRQ(ierr);
621         ierr = KSPGetType(origksp,&ksp_type);CHKERRQ(ierr);
622         ierr = KSPSetType(schurksp,ksp_type);CHKERRQ(ierr);
623         ierr = KSPGetPC(schurksp,&schurpc);CHKERRQ(ierr);
624         ierr = KSPGetPC(origksp,&origpc);CHKERRQ(ierr);
625         ierr = PetscObjectTypeCompare((PetscObject)origpc,PCNONE,&ispcnone);CHKERRQ(ierr);
626         if (!ispcnone) {
627           PCType pc_type;
628           ierr = PCGetType(origpc,&pc_type);CHKERRQ(ierr);
629           ierr = PCSetType(schurpc,pc_type);CHKERRQ(ierr);
630         } else {
631           ierr = PCSetType(schurpc,PCLU);CHKERRQ(ierr);
632         }
633         ierr = ISGetSize(is_I,&n_internal);CHKERRQ(ierr);
634         if (n_internal) { /* UMFPACK gives error with 0 sized problems */
635           MatSolverPackage solver=NULL;
636           ierr = PCFactorGetMatSolverPackage(origpc,(const MatSolverPackage*)&solver);CHKERRQ(ierr);
637           if (solver) {
638             ierr = PCFactorSetMatSolverPackage(schurpc,solver);CHKERRQ(ierr);
639           }
640         }
641         ierr = KSPSetUp(schurksp);CHKERRQ(ierr);
642       }
643       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
644       ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&S_Ej_expl);CHKERRQ(ierr);
645       ierr = PCBDDCComputeExplicitSchur(S_Ej,sub_schurs->is_hermitian,MAT_REUSE_MATRIX,&S_Ej_expl);CHKERRQ(ierr);
646       ierr = PetscObjectTypeCompare((PetscObject)S_Ej_expl,MATSEQDENSE,&Sdense);CHKERRQ(ierr);
647       if (Sdense) {
648         for (j=0;j<subset_size;j++) {
649           dummy_idx[j]=local_size+j;
650         }
651         ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
652       } else {
653         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not yet implemented for sparse matrices");
654       }
655       ierr = MatDestroy(&S_Ej);CHKERRQ(ierr);
656       ierr = MatDestroy(&S_Ej_expl);CHKERRQ(ierr);
657       local_size += subset_size;
658     }
659     ierr = PetscFree2(dummy_idx,work);CHKERRQ(ierr);
660     /* free */
661     ierr = ISDestroy(&is_I);CHKERRQ(ierr);
662     ierr = MatDestroy(&AE_II);CHKERRQ(ierr);
663     ierr = PetscFree(all_local_idx_N);CHKERRQ(ierr);
664   } else {
665     Mat         A;
666     IS          is_A_all;
667     PetscScalar *work,*S_data,*schur_factor;
668     PetscInt    n,n_I,*dummy_idx,size_schur,size_active_schur,cum,cum2;
669     PetscBool   economic,mumps_S,S_lower_triangular = PETSC_FALSE,factor_workaround;
670 
671     /* get sizes */
672     n_I = 0;
673     if (is_I_layer) {
674       ierr = ISGetLocalSize(is_I_layer,&n_I);CHKERRQ(ierr);
675     }
676     economic = PETSC_FALSE;
677     ierr = ISGetLocalSize(sub_schurs->is_I,&cum);CHKERRQ(ierr);
678     if (cum != n_I) economic = PETSC_TRUE;
679     ierr = MatGetLocalSize(sub_schurs->A,&n,NULL);CHKERRQ(ierr);
680 
681     /* size active schurs does not count any dirichlet or vertex dof on the interface */
682     size_active_schur = local_size;
683     cum = n_I+size_active_schur;
684     if (sub_schurs->is_dir) {
685       const PetscInt* idxs;
686       PetscInt        n_dir;
687 
688       ierr = ISGetLocalSize(sub_schurs->is_dir,&n_dir);CHKERRQ(ierr);
689       ierr = ISGetIndices(sub_schurs->is_dir,&idxs);CHKERRQ(ierr);
690       ierr = PetscMemcpy(all_local_idx_N+cum,idxs,n_dir*sizeof(PetscInt));CHKERRQ(ierr);
691       ierr = ISRestoreIndices(sub_schurs->is_dir,&idxs);CHKERRQ(ierr);
692       cum += n_dir;
693     }
694     factor_workaround = PETSC_FALSE;
695     /* include the primal vertices in the Schur complement */
696     if (exact_schur && sub_schurs->is_vertices && compute_Stilda) {
697       PetscInt n_v;
698 
699       ierr = ISGetLocalSize(sub_schurs->is_vertices,&n_v);CHKERRQ(ierr);
700       if (n_v) {
701         const PetscInt* idxs;
702 
703         ierr = ISGetIndices(sub_schurs->is_vertices,&idxs);CHKERRQ(ierr);
704         ierr = PetscMemcpy(all_local_idx_N+cum,idxs,n_v*sizeof(PetscInt));CHKERRQ(ierr);
705         ierr = ISRestoreIndices(sub_schurs->is_vertices,&idxs);CHKERRQ(ierr);
706         cum += n_v;
707         factor_workaround = PETSC_TRUE;
708       }
709     }
710     size_schur = cum - n_I;
711     ierr = ISCreateGeneral(PETSC_COMM_SELF,cum,all_local_idx_N,PETSC_OWN_POINTER,&is_A_all);CHKERRQ(ierr);
712     /* get working mat (TODO: factorize without actually permuting)  */
713     if (cum == n) {
714       ierr = ISSetPermutation(is_A_all);CHKERRQ(ierr);
715       ierr = MatPermute(sub_schurs->A,is_A_all,is_A_all,&A);CHKERRQ(ierr);
716     } else {
717       ierr = MatGetSubMatrix(sub_schurs->A,is_A_all,is_A_all,MAT_INITIAL_MATRIX,&A);CHKERRQ(ierr);
718     }
719     ierr = MatSetOptionsPrefix(A,"sub_schurs_");CHKERRQ(ierr);
720     ierr = MatSetOption(A,MAT_SYMMETRIC,sub_schurs->is_hermitian);CHKERRQ(ierr);
721     ierr = MatSetOption(A,MAT_HERMITIAN,sub_schurs->is_hermitian);CHKERRQ(ierr);
722     ierr = MatSetOption(A,MAT_SPD,sub_schurs->is_posdef);CHKERRQ(ierr);
723 
724     /* when using the benign subspace trick, the local Schur complements are SPD */
725     if (benign_trick) sub_schurs->is_posdef = PETSC_TRUE;
726 
727     if (n_I) { /* TODO add ordering from options */
728       IS is_schur;
729 
730       if (sub_schurs->is_hermitian && sub_schurs->is_posdef) {
731         ierr = MatGetFactor(A,MATSOLVERMUMPS,MAT_FACTOR_CHOLESKY,&F);CHKERRQ(ierr);
732       } else {
733         ierr = MatGetFactor(A,MATSOLVERMUMPS,MAT_FACTOR_LU,&F);CHKERRQ(ierr);
734       }
735       /* subsets ordered last */
736       ierr = ISCreateStride(PETSC_COMM_SELF,size_schur,n_I,1,&is_schur);CHKERRQ(ierr);
737       ierr = MatFactorSetSchurIS(F,is_schur);CHKERRQ(ierr);
738       ierr = ISDestroy(&is_schur);CHKERRQ(ierr);
739 
740       /* factorization step */
741       if (sub_schurs->is_hermitian && sub_schurs->is_posdef) {
742         ierr = MatCholeskyFactorSymbolic(F,A,NULL,NULL);CHKERRQ(ierr);
743 #if defined(PETSC_HAVE_MUMPS) /* be sure that icntl 19 is not set by command line */
744         ierr = MatMumpsSetIcntl(F,19,2);CHKERRQ(ierr);
745         if (sub_schurs->is_posdef) {
746           ierr = MatMumpsSchurComplementSetSym(F,1);CHKERRQ(ierr);
747         } else {
748           ierr = MatMumpsSchurComplementSetSym(F,2);CHKERRQ(ierr);
749         }
750 #endif
751         ierr = MatCholeskyFactorNumeric(F,A,NULL);CHKERRQ(ierr);
752         S_lower_triangular = PETSC_TRUE;
753       } else {
754         ierr = MatLUFactorSymbolic(F,A,NULL,NULL,NULL);CHKERRQ(ierr);
755 #if defined(PETSC_HAVE_MUMPS) /* be sure that icntl 19 is not set by command line */
756         ierr = MatMumpsSetIcntl(F,19,3);CHKERRQ(ierr);
757 #endif
758         ierr = MatLUFactorNumeric(F,A,NULL);CHKERRQ(ierr);
759       }
760 
761       /* get explicit Schur Complement computed during numeric factorization */
762       ierr = MatFactorGetSchurComplement(F,&S_all);CHKERRQ(ierr);
763 
764       /* we can reuse the solvers if we are not using the economic version */
765       reuse_solvers = (PetscBool)(reuse_solvers && !economic);
766       factor_workaround = (PetscBool)(reuse_solvers && factor_workaround);
767       mumps_S = PETSC_TRUE;
768     } else { /* we can't use MUMPS when size_schur == size_of_the_problem */
769       ierr = MatConvert(A,MATSEQDENSE,MAT_INITIAL_MATRIX,&S_all);CHKERRQ(ierr);
770       reuse_solvers = PETSC_FALSE; /* TODO: why we can't reuse the solvers here? */
771       mumps_S = PETSC_FALSE;
772     }
773 
774     if (reuse_solvers) {
775       Mat              A_II;
776       Vec              vec1_B;
777       PCBDDCReuseMumps msolv_ctx;
778 
779       if (sub_schurs->reuse_mumps) {
780         ierr = PCBDDCReuseMumpsReset(sub_schurs->reuse_mumps);CHKERRQ(ierr);
781       } else {
782         ierr = PetscNew(&sub_schurs->reuse_mumps);CHKERRQ(ierr);
783       }
784       msolv_ctx = sub_schurs->reuse_mumps;
785       ierr = MatSchurComplementGetSubMatrices(sub_schurs->S,&A_II,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
786       ierr = PetscObjectReference((PetscObject)F);CHKERRQ(ierr);
787       msolv_ctx->F = F;
788       ierr = MatCreateVecs(F,&msolv_ctx->sol,NULL);CHKERRQ(ierr);
789       /* currently PETSc has no support for MatSolve(F,x,x), so cheat and let rhs and sol share the same memory */
790       {
791         PetscScalar *array;
792         PetscInt    n;
793 
794         ierr = VecGetLocalSize(msolv_ctx->sol,&n);CHKERRQ(ierr);
795         ierr = VecGetArray(msolv_ctx->sol,&array);CHKERRQ(ierr);
796         ierr = VecCreateSeqWithArray(PetscObjectComm((PetscObject)msolv_ctx->sol),1,n,array,&msolv_ctx->rhs);CHKERRQ(ierr);
797         ierr = VecRestoreArray(msolv_ctx->sol,&array);CHKERRQ(ierr);
798       }
799       msolv_ctx->has_vertices = factor_workaround;
800 
801       /* interior solver */
802       ierr = PCCreate(PetscObjectComm((PetscObject)A_II),&msolv_ctx->interior_solver);CHKERRQ(ierr);
803       ierr = PCSetOperators(msolv_ctx->interior_solver,A_II,A_II);CHKERRQ(ierr);
804       ierr = PCSetType(msolv_ctx->interior_solver,PCSHELL);CHKERRQ(ierr);
805       ierr = PCShellSetContext(msolv_ctx->interior_solver,msolv_ctx);CHKERRQ(ierr);
806       ierr = PCShellSetApply(msolv_ctx->interior_solver,PCBDDCMumpsInteriorSolve);CHKERRQ(ierr);
807       ierr = PCShellSetApplyTranspose(msolv_ctx->interior_solver,PCBDDCMumpsInteriorSolveTranspose);CHKERRQ(ierr);
808 
809       /* correction solver */
810       ierr = PCCreate(PetscObjectComm((PetscObject)A),&msolv_ctx->correction_solver);CHKERRQ(ierr);
811       ierr = PCSetOperators(msolv_ctx->correction_solver,A,A);CHKERRQ(ierr);
812       ierr = PCSetType(msolv_ctx->correction_solver,PCSHELL);CHKERRQ(ierr);
813       ierr = PCShellSetContext(msolv_ctx->correction_solver,msolv_ctx);CHKERRQ(ierr);
814       ierr = PCShellSetApply(msolv_ctx->correction_solver,PCBDDCMumpsCorrectionSolve);CHKERRQ(ierr);
815       ierr = PCShellSetApplyTranspose(msolv_ctx->correction_solver,PCBDDCMumpsCorrectionSolveTranspose);CHKERRQ(ierr);
816 
817       /* scatter and vecs for Schur complement solver */
818       ierr = MatCreateVecs(S_all,&msolv_ctx->sol_B,&msolv_ctx->rhs_B);CHKERRQ(ierr);
819       ierr = MatCreateVecs(sub_schurs->S,&vec1_B,NULL);CHKERRQ(ierr);
820       if (!factor_workaround) {
821         ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,is_A_all,&msolv_ctx->is_B);CHKERRQ(ierr);
822         ierr = VecScatterCreate(vec1_B,msolv_ctx->is_B,msolv_ctx->sol_B,NULL,&msolv_ctx->correction_scatter_B);CHKERRQ(ierr);
823         ierr = PetscObjectReference((PetscObject)is_A_all);CHKERRQ(ierr);
824         msolv_ctx->is_R = is_A_all;
825       } else {
826         IS              is_B_all;
827         const PetscInt* idxs;
828         PetscInt        dual,n_v,n;
829 
830         ierr = ISGetLocalSize(sub_schurs->is_vertices,&n_v);CHKERRQ(ierr);
831         dual = size_schur - n_v;
832         ierr = ISGetLocalSize(is_A_all,&n);CHKERRQ(ierr);
833         ierr = ISGetIndices(is_A_all,&idxs);CHKERRQ(ierr);
834         ierr = ISCreateGeneral(PetscObjectComm((PetscObject)is_A_all),dual,idxs+n_I,PETSC_COPY_VALUES,&is_B_all);CHKERRQ(ierr);
835         ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,is_B_all,&msolv_ctx->is_B);CHKERRQ(ierr);
836         ierr = ISDestroy(&is_B_all);CHKERRQ(ierr);
837         ierr = ISCreateStride(PetscObjectComm((PetscObject)is_A_all),dual,0,1,&is_B_all);CHKERRQ(ierr);
838         ierr = VecScatterCreate(vec1_B,msolv_ctx->is_B,msolv_ctx->sol_B,is_B_all,&msolv_ctx->correction_scatter_B);CHKERRQ(ierr);
839         ierr = ISDestroy(&is_B_all);CHKERRQ(ierr);
840         ierr = ISCreateGeneral(PetscObjectComm((PetscObject)is_A_all),n-n_v,idxs,PETSC_COPY_VALUES,&msolv_ctx->is_R);CHKERRQ(ierr);
841         ierr = ISRestoreIndices(is_A_all,&idxs);CHKERRQ(ierr);
842       }
843       ierr = VecDestroy(&vec1_B);CHKERRQ(ierr);
844     }
845     ierr = MatDestroy(&A);CHKERRQ(ierr);
846     ierr = ISDestroy(&is_A_all);CHKERRQ(ierr);
847 
848     /* Work arrays */
849     ierr = PetscMalloc2(max_subset_size,&dummy_idx,max_subset_size*max_subset_size,&work);CHKERRQ(ierr);
850 
851     /* matrices for adaptive selection */
852     if (compute_Stilda) {
853       if (!sub_schurs->sum_S_Ej_tilda_all) {
854         ierr = MatCreate(PETSC_COMM_SELF,&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
855         ierr = MatSetSizes(sub_schurs->sum_S_Ej_tilda_all,PETSC_DECIDE,PETSC_DECIDE,size_active_schur,size_active_schur);CHKERRQ(ierr);
856         ierr = MatSetType(sub_schurs->sum_S_Ej_tilda_all,MATAIJ);CHKERRQ(ierr);
857         ierr = MatSeqAIJSetPreallocation(sub_schurs->sum_S_Ej_tilda_all,0,nnz);CHKERRQ(ierr);
858       }
859       if (!sub_schurs->sum_S_Ej_inv_all) {
860         ierr = MatCreate(PETSC_COMM_SELF,&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
861         ierr = MatSetSizes(sub_schurs->sum_S_Ej_inv_all,PETSC_DECIDE,PETSC_DECIDE,size_active_schur,size_active_schur);CHKERRQ(ierr);
862         ierr = MatSetType(sub_schurs->sum_S_Ej_inv_all,MATAIJ);CHKERRQ(ierr);
863         ierr = MatSeqAIJSetPreallocation(sub_schurs->sum_S_Ej_inv_all,0,nnz);CHKERRQ(ierr);
864       }
865     }
866 
867     /* S_Ej_all */
868     cum = cum2 = 0;
869     ierr = MatDenseGetArray(S_all,&S_data);CHKERRQ(ierr);
870     for (i=0;i<sub_schurs->n_subs;i++) {
871       PetscInt j;
872 
873       /* get S_E */
874       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
875       if (S_lower_triangular) { /* I need to expand the upper triangular data (column oriented) */
876         PetscInt k;
877         for (k=0;k<subset_size;k++) {
878           for (j=k;j<subset_size;j++) {
879             work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
880             work[j*subset_size+k] = PetscConj(S_data[cum2+k*size_schur+j]);
881           }
882         }
883       } else { /* just copy to workspace */
884         PetscInt k;
885         for (k=0;k<subset_size;k++) {
886           for (j=0;j<subset_size;j++) {
887             work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
888           }
889         }
890       }
891       /* insert S_E values */
892       for (j=0;j<subset_size;j++) dummy_idx[j] = cum+j;
893       ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
894 
895       /* if adaptivity is requested, invert S_E blocks */
896       if (compute_Stilda) {
897         ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr);
898         ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
899         if (!use_getr) { /* TODO add sytrf/i for symmetric non hermitian */
900           PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,work,&B_N,&B_ierr));
901           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
902           PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,work,&B_N,&B_ierr));
903           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
904         } else {
905           PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,work,&B_N,pivots,&B_ierr));
906           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
907           PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,work,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
908           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
909         }
910         ierr = PetscFPTrapPop();CHKERRQ(ierr);
911         ierr = MatSetValues(sub_schurs->sum_S_Ej_inv_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
912       }
913       cum += subset_size;
914       cum2 += subset_size*(size_schur + 1);
915     }
916     ierr = MatDenseRestoreArray(S_all,&S_data);CHKERRQ(ierr);
917 
918     if (mumps_S) {
919       ierr = MatFactorRestoreSchurComplement(F,&S_all);CHKERRQ(ierr);
920     }
921 
922     schur_factor = NULL;
923     if (compute_Stilda && size_active_schur) {
924 
925       if (sub_schurs->n_subs == 1 && size_schur == size_active_schur) { /* we already computed the inverse */
926         PetscInt j;
927         for (j=0;j<size_schur;j++) dummy_idx[j] = j;
928         ierr = MatSetValues(sub_schurs->sum_S_Ej_tilda_all,size_schur,dummy_idx,size_schur,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
929       } else {
930         Mat S_all_inv=NULL;
931         if (mumps_S) { /* use MatFactor calls to invert S */
932             /* let MatFactor factorize the Schur complement */
933           ierr = MatFactorFactorizeSchurComplement(F);CHKERRQ(ierr);
934           /* for adaptive selection we need S^-1; for solver reusage we need S_\Delta\Delta^-1.
935              The latter is not the principal subminor for S^-1. However, the factors can be reused since S_\Delta\Delta is the leading principal submatrix of S */
936           if (factor_workaround) {
937             PetscScalar *data;
938             PetscInt nd = 0;
939 
940             if (!sub_schurs->is_posdef) {
941               SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor update not yet implemented for non SPD matrices");
942             }
943             ierr = MatFactorGetSchurComplement(F,&S_all_inv);CHKERRQ(ierr);
944             ierr = MatDenseGetArray(S_all_inv,&data);CHKERRQ(ierr);
945             if (sub_schurs->is_dir) { /* dirichlet dofs could have different scalings */
946               ierr = ISGetLocalSize(sub_schurs->is_dir,&nd);CHKERRQ(ierr);
947             }
948             ierr = PetscMalloc1((size_active_schur*(size_active_schur+1))/2+nd,&schur_factor);CHKERRQ(ierr);
949             cum = 0;
950             for (i=0;i<size_active_schur;i++) {
951               ierr = PetscMemcpy(schur_factor+cum,data+i*(size_schur+1),(size_active_schur-i)*sizeof(PetscScalar));CHKERRQ(ierr);
952               cum += size_active_schur-i;
953             }
954             for (i=0;i<nd;i++) schur_factor[cum+i] = data[(i+size_active_schur)*(size_schur+1)];
955             /* invert without calling MatFactorInvertSchurComplement, since we are hacking */
956             ierr = PetscBLASIntCast(size_schur,&B_N);CHKERRQ(ierr);
957             ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
958             PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,data,&B_N,&B_ierr));
959             ierr = PetscFPTrapPop();CHKERRQ(ierr);
960             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
961             ierr = MatDenseRestoreArray(S_all_inv,&data);CHKERRQ(ierr);
962           } else {
963             ierr = MatFactorInvertSchurComplement(F);CHKERRQ(ierr);
964             ierr = MatFactorGetSchurComplement(F,&S_all_inv);CHKERRQ(ierr);
965           }
966         } else { /* we need to invert explicitly since we are not using MUMPS for S */
967           ierr = PetscObjectReference((PetscObject)S_all);CHKERRQ(ierr);
968           S_all_inv = S_all;
969           ierr = MatDenseGetArray(S_all_inv,&S_data);CHKERRQ(ierr);
970           ierr = PetscBLASIntCast(size_schur,&B_N);CHKERRQ(ierr);
971           ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
972           if (!use_getr) { /* TODO add sytrf/i for symmetric non hermitian */
973             PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,S_data,&B_N,&B_ierr));
974             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
975             PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,S_data,&B_N,&B_ierr));
976             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
977           } else {
978             PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,S_data,&B_N,pivots,&B_ierr));
979             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
980             PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,S_data,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
981             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
982           }
983           ierr = PetscFPTrapPop();CHKERRQ(ierr);
984           ierr = MatDenseRestoreArray(S_all_inv,&S_data);CHKERRQ(ierr);
985         }
986         /* S_Ej_tilda_all */
987         cum = cum2 = 0;
988         ierr = MatDenseGetArray(S_all_inv,&S_data);CHKERRQ(ierr);
989         for (i=0;i<sub_schurs->n_subs;i++) {
990           PetscInt j;
991 
992           ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
993           /* get (St^-1)_E */
994           /* Here I don't need to expand to upper triangular since St^-1
995              will be properly accessed later during adaptive selection */
996           if (S_lower_triangular) {
997             PetscInt k;
998             for (k=0;k<subset_size;k++) {
999               for (j=k;j<subset_size;j++) {
1000                 work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
1001               }
1002             }
1003           } else { /* copy to workspace */
1004             PetscInt k;
1005             for (k=0;k<subset_size;k++) {
1006               for (j=0;j<subset_size;j++) {
1007                 work[k*subset_size+j] = S_data[cum2+k*size_schur+j];
1008               }
1009             }
1010           }
1011           for (j=0;j<subset_size;j++) dummy_idx[j] = cum+j;
1012           ierr = MatSetValues(sub_schurs->sum_S_Ej_tilda_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
1013           cum += subset_size;
1014           cum2 += subset_size*(size_schur + 1);
1015         }
1016         ierr = MatDenseRestoreArray(S_all_inv,&S_data);CHKERRQ(ierr);
1017         if (mumps_S) {
1018           ierr = MatFactorRestoreSchurComplement(F,&S_all_inv);CHKERRQ(ierr);
1019         }
1020         ierr = MatDestroy(&S_all_inv);CHKERRQ(ierr);
1021       }
1022 
1023       /* move back factors to Schur data into F */
1024       if (factor_workaround) {
1025         Mat S_tmp;
1026         PetscInt nv,nd=0;
1027 
1028         if (!mumps_S) {
1029           SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"This should not happen");
1030         }
1031         ierr = ISGetLocalSize(sub_schurs->is_vertices,&nv);CHKERRQ(ierr);
1032         ierr = MatFactorGetSchurComplement(F,&S_tmp);CHKERRQ(ierr);
1033         if (sub_schurs->is_posdef) {
1034           PetscScalar *data;
1035 
1036           ierr = MatZeroEntries(S_tmp);CHKERRQ(ierr);
1037           ierr = MatDenseGetArray(S_tmp,&data);CHKERRQ(ierr);
1038 
1039           if (S_lower_triangular) {
1040             cum = 0;
1041             for (i=0;i<size_active_schur;i++) {
1042               ierr = PetscMemcpy(data+i*(size_schur+1),schur_factor+cum,(size_active_schur-i)*sizeof(PetscScalar));CHKERRQ(ierr);
1043               cum += size_active_schur-i;
1044 	    }
1045           } else {
1046             ierr = PetscMemcpy(data,schur_factor,size_schur*size_schur*sizeof(PetscScalar));CHKERRQ(ierr);
1047           }
1048           if (sub_schurs->is_dir) {
1049             ierr = ISGetLocalSize(sub_schurs->is_dir,&nd);CHKERRQ(ierr);
1050 	    for (i=0;i<nd;i++) {
1051 	      data[(i+size_active_schur)*(size_schur+1)] = schur_factor[cum+i];
1052 	    }
1053           }
1054           /* workaround: since I cannot modify the matrices used inside the solvers for the forward and backward substitutions,
1055              set the diagonal entry of the Schur factor to a very large value */
1056           for (i=size_active_schur+nd;i<size_schur;i++) {
1057             data[i*(size_schur+1)] = PETSC_MAX_REAL;
1058           }
1059           ierr = MatDenseRestoreArray(S_tmp,&data);CHKERRQ(ierr);
1060         } else {
1061           SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor update not yet implemented for non SPD matrices");
1062         }
1063         ierr = MatFactorRestoreSchurComplement(F,&S_tmp);CHKERRQ(ierr);
1064       }
1065     }
1066     ierr = PetscFree2(dummy_idx,work);CHKERRQ(ierr);
1067     ierr = PetscFree(schur_factor);CHKERRQ(ierr);
1068   }
1069   ierr = PetscFree(nnz);CHKERRQ(ierr);
1070   ierr = MatDestroy(&F);CHKERRQ(ierr);
1071   ierr = ISDestroy(&is_I_layer);CHKERRQ(ierr);
1072   ierr = MatDestroy(&S_all);CHKERRQ(ierr);
1073   ierr = MatDestroy(&A_BB);CHKERRQ(ierr);
1074   ierr = MatDestroy(&A_IB);CHKERRQ(ierr);
1075   ierr = MatDestroy(&A_BI);CHKERRQ(ierr);
1076   ierr = MatAssemblyBegin(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1077   ierr = MatAssemblyEnd(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1078   if (compute_Stilda) {
1079     ierr = MatAssemblyBegin(sub_schurs->sum_S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1080     ierr = MatAssemblyEnd(sub_schurs->sum_S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1081     ierr = MatAssemblyBegin(sub_schurs->sum_S_Ej_inv_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1082     ierr = MatAssemblyEnd(sub_schurs->sum_S_Ej_inv_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1083   }
1084 
1085   /* Global matrix of all assembled Schur on subsets */
1086   ierr = MatISSetLocalMat(work_mat,sub_schurs->S_Ej_all);CHKERRQ(ierr);
1087   ierr = MatISSetMPIXAIJPreallocation_Private(work_mat,global_schur_subsets,PETSC_TRUE);CHKERRQ(ierr);
1088   ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
1089 
1090   /* Get local part of (\sum_j S_Ej) */
1091   if (!sub_schurs->sum_S_Ej_all) {
1092     ierr = MatGetSubMatrices(global_schur_subsets,1,&all_subsets_n,&all_subsets_n,MAT_INITIAL_MATRIX,&submats);CHKERRQ(ierr);
1093     sub_schurs->sum_S_Ej_all = submats[0];
1094   } else {
1095     ierr = PetscMalloc1(1,&submats);CHKERRQ(ierr);
1096     submats[0] = sub_schurs->sum_S_Ej_all;
1097     ierr = MatGetSubMatrices(global_schur_subsets,1,&all_subsets_n,&all_subsets_n,MAT_REUSE_MATRIX,&submats);CHKERRQ(ierr);
1098   }
1099 
1100   /* Compute explicitly (\sum_j S_Ej)^-1 (faster scaling during PCApply, needs extra work when doing setup) */
1101   if (faster_deluxe) {
1102     Mat         tmpmat;
1103     PetscScalar *array;
1104     PetscInt    cum;
1105 
1106     ierr = MatSeqAIJGetArray(sub_schurs->sum_S_Ej_all,&array);CHKERRQ(ierr);
1107     cum = 0;
1108     for (i=0;i<sub_schurs->n_subs;i++) {
1109       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
1110       ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr);
1111       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
1112       if (!use_getr) {
1113         PetscInt j,k;
1114 
1115         PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,array+cum,&B_N,&B_ierr));
1116         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
1117         PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,array+cum,&B_N,&B_ierr));
1118         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
1119         for (j=0;j<B_N;j++) {
1120           for (k=j+1;k<B_N;k++) {
1121             array[k*B_N+j+cum] = array[j*B_N+k+cum];
1122           }
1123         }
1124       } else {
1125         PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,array+cum,&B_N,pivots,&B_ierr));
1126         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
1127         PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,array+cum,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
1128         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
1129       }
1130       ierr = PetscFPTrapPop();CHKERRQ(ierr);
1131       cum += subset_size*subset_size;
1132     }
1133     ierr = MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_all,&array);CHKERRQ(ierr);
1134     ierr = MatMatMult(sub_schurs->S_Ej_all,sub_schurs->sum_S_Ej_all,MAT_INITIAL_MATRIX,1.0,&tmpmat);CHKERRQ(ierr);
1135     ierr = MatDestroy(&sub_schurs->S_Ej_all);CHKERRQ(ierr);
1136     ierr = MatDestroy(&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
1137     sub_schurs->S_Ej_all = tmpmat;
1138   }
1139 
1140   /* Get local part of (\sum_j S^-1_Ej) (\sum_j St^-1_Ej) */
1141   if (compute_Stilda) {
1142     ierr = MatISSetLocalMat(work_mat,sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
1143     ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
1144     submats[0] = sub_schurs->sum_S_Ej_tilda_all;
1145     ierr = MatGetSubMatrices(global_schur_subsets,1,&all_subsets_n,&all_subsets_n,MAT_REUSE_MATRIX,&submats);CHKERRQ(ierr);
1146     ierr = MatISSetLocalMat(work_mat,sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
1147     ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
1148     submats[0] = sub_schurs->sum_S_Ej_inv_all;
1149     ierr = MatGetSubMatrices(global_schur_subsets,1,&all_subsets_n,&all_subsets_n,MAT_REUSE_MATRIX,&submats);CHKERRQ(ierr);
1150   }
1151 
1152   /* free workspace */
1153   ierr = PetscFree(submats);CHKERRQ(ierr);
1154   ierr = PetscFree2(Bwork,pivots);CHKERRQ(ierr);
1155   ierr = MatDestroy(&global_schur_subsets);CHKERRQ(ierr);
1156   ierr = MatDestroy(&work_mat);CHKERRQ(ierr);
1157   ierr = ISDestroy(&all_subsets_n);CHKERRQ(ierr);
1158   ierr = PetscCommDestroy(&comm_n);CHKERRQ(ierr);
1159   PetscFunctionReturn(0);
1160 }
1161 
1162 #undef __FUNCT__
1163 #define __FUNCT__ "PCBDDCSubSchursInit"
1164 PetscErrorCode PCBDDCSubSchursInit(PCBDDCSubSchurs sub_schurs, IS is_I, IS is_B, PCBDDCGraph graph, ISLocalToGlobalMapping BtoNmap)
1165 {
1166   IS              *faces,*edges,*all_cc,vertices;
1167   PetscInt        i,n_faces,n_edges,n_all_cc;
1168   PetscBool       is_sorted;
1169   PetscErrorCode  ierr;
1170 
1171   PetscFunctionBegin;
1172   ierr = ISSorted(is_I,&is_sorted);CHKERRQ(ierr);
1173   if (!is_sorted) {
1174     SETERRQ(PetscObjectComm((PetscObject)is_I),PETSC_ERR_PLIB,"IS for I dofs should be shorted");
1175   }
1176   ierr = ISSorted(is_B,&is_sorted);CHKERRQ(ierr);
1177   if (!is_sorted) {
1178     SETERRQ(PetscObjectComm((PetscObject)is_B),PETSC_ERR_PLIB,"IS for B dofs should be shorted");
1179   }
1180 
1181   /* reset any previous data */
1182   ierr = PCBDDCSubSchursReset(sub_schurs);CHKERRQ(ierr);
1183 
1184   /* get index sets for faces and edges (already sorted by global ordering) */
1185   ierr = PCBDDCGraphGetCandidatesIS(graph,&n_faces,&faces,&n_edges,&edges,&vertices);CHKERRQ(ierr);
1186   n_all_cc = n_faces+n_edges;
1187   ierr = PetscBTCreate(n_all_cc,&sub_schurs->is_edge);CHKERRQ(ierr);
1188   ierr = PetscMalloc1(n_all_cc,&all_cc);CHKERRQ(ierr);
1189   for (i=0;i<n_faces;i++) {
1190     all_cc[i] = faces[i];
1191   }
1192   for (i=0;i<n_edges;i++) {
1193     all_cc[n_faces+i] = edges[i];
1194     ierr = PetscBTSet(sub_schurs->is_edge,n_faces+i);CHKERRQ(ierr);
1195   }
1196   ierr = PetscFree(faces);CHKERRQ(ierr);
1197   ierr = PetscFree(edges);CHKERRQ(ierr);
1198   sub_schurs->is_dir = NULL;
1199   ierr = PCBDDCGraphGetDirichletDofsB(graph,&sub_schurs->is_dir);CHKERRQ(ierr);
1200 
1201   /* Determine if MUMPS can be used */
1202   sub_schurs->use_mumps = PETSC_FALSE;
1203 #if defined(PETSC_HAVE_MUMPS)
1204   sub_schurs->use_mumps = PETSC_TRUE;
1205 #endif
1206 
1207   ierr = PetscObjectReference((PetscObject)is_I);CHKERRQ(ierr);
1208   sub_schurs->is_I = is_I;
1209   ierr = PetscObjectReference((PetscObject)is_B);CHKERRQ(ierr);
1210   sub_schurs->is_B = is_B;
1211   ierr = PetscObjectReference((PetscObject)graph->l2gmap);CHKERRQ(ierr);
1212   sub_schurs->l2gmap = graph->l2gmap;
1213   ierr = PetscObjectReference((PetscObject)BtoNmap);CHKERRQ(ierr);
1214   sub_schurs->BtoNmap = BtoNmap;
1215   sub_schurs->n_subs = n_all_cc;
1216   sub_schurs->is_subs = all_cc;
1217   if (!sub_schurs->use_mumps) { /* sort by local ordering mumps is not present */
1218     for (i=0;i<sub_schurs->n_subs;i++) {
1219       ierr = ISSort(sub_schurs->is_subs[i]);CHKERRQ(ierr);
1220     }
1221   }
1222   sub_schurs->is_vertices = vertices;
1223   sub_schurs->S_Ej_all = NULL;
1224   sub_schurs->sum_S_Ej_all = NULL;
1225   sub_schurs->sum_S_Ej_inv_all = NULL;
1226   sub_schurs->sum_S_Ej_tilda_all = NULL;
1227   sub_schurs->is_Ej_all = NULL;
1228   PetscFunctionReturn(0);
1229 }
1230 
1231 #undef __FUNCT__
1232 #define __FUNCT__ "PCBDDCSubSchursCreate"
1233 PetscErrorCode PCBDDCSubSchursCreate(PCBDDCSubSchurs *sub_schurs)
1234 {
1235   PCBDDCSubSchurs schurs_ctx;
1236   PetscErrorCode  ierr;
1237 
1238   PetscFunctionBegin;
1239   ierr = PetscNew(&schurs_ctx);CHKERRQ(ierr);
1240   schurs_ctx->n_subs = 0;
1241   *sub_schurs = schurs_ctx;
1242   PetscFunctionReturn(0);
1243 }
1244 
1245 #undef __FUNCT__
1246 #define __FUNCT__ "PCBDDCSubSchursDestroy"
1247 PetscErrorCode PCBDDCSubSchursDestroy(PCBDDCSubSchurs *sub_schurs)
1248 {
1249   PetscErrorCode ierr;
1250 
1251   PetscFunctionBegin;
1252   ierr = PCBDDCSubSchursReset(*sub_schurs);CHKERRQ(ierr);
1253   ierr = PetscFree(*sub_schurs);CHKERRQ(ierr);
1254   PetscFunctionReturn(0);
1255 }
1256 
1257 #undef __FUNCT__
1258 #define __FUNCT__ "PCBDDCSubSchursReset"
1259 PetscErrorCode PCBDDCSubSchursReset(PCBDDCSubSchurs sub_schurs)
1260 {
1261   PetscInt       i;
1262   PetscErrorCode ierr;
1263 
1264   PetscFunctionBegin;
1265   ierr = MatDestroy(&sub_schurs->A);CHKERRQ(ierr);
1266   ierr = MatDestroy(&sub_schurs->S);CHKERRQ(ierr);
1267   ierr = ISDestroy(&sub_schurs->is_I);CHKERRQ(ierr);
1268   ierr = ISDestroy(&sub_schurs->is_B);CHKERRQ(ierr);
1269   ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->l2gmap);CHKERRQ(ierr);
1270   ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->BtoNmap);CHKERRQ(ierr);
1271   ierr = MatDestroy(&sub_schurs->S_Ej_all);CHKERRQ(ierr);
1272   ierr = MatDestroy(&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
1273   ierr = MatDestroy(&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
1274   ierr = MatDestroy(&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
1275   ierr = ISDestroy(&sub_schurs->is_Ej_all);CHKERRQ(ierr);
1276   ierr = ISDestroy(&sub_schurs->is_vertices);CHKERRQ(ierr);
1277   ierr = ISDestroy(&sub_schurs->is_dir);CHKERRQ(ierr);
1278   ierr = PetscBTDestroy(&sub_schurs->is_edge);CHKERRQ(ierr);
1279   for (i=0;i<sub_schurs->n_subs;i++) {
1280     ierr = ISDestroy(&sub_schurs->is_subs[i]);CHKERRQ(ierr);
1281   }
1282   if (sub_schurs->n_subs) {
1283     ierr = PetscFree(sub_schurs->is_subs);CHKERRQ(ierr);
1284   }
1285   if (sub_schurs->reuse_mumps) {
1286     ierr = PCBDDCReuseMumpsReset(sub_schurs->reuse_mumps);CHKERRQ(ierr);
1287   }
1288   ierr = PetscFree(sub_schurs->reuse_mumps);CHKERRQ(ierr);
1289   sub_schurs->n_subs = 0;
1290   PetscFunctionReturn(0);
1291 }
1292 
1293 #undef __FUNCT__
1294 #define __FUNCT__ "PCBDDCAdjGetNextLayer_Private"
1295 PETSC_STATIC_INLINE PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt* queue_tip,PetscInt n_prev,PetscBT touched,PetscInt* xadj,PetscInt* adjncy,PetscInt* n_added)
1296 {
1297   PetscInt       i,j,n;
1298   PetscErrorCode ierr;
1299 
1300   PetscFunctionBegin;
1301   n = 0;
1302   for (i=-n_prev;i<0;i++) {
1303     PetscInt start_dof = queue_tip[i];
1304     for (j=xadj[start_dof];j<xadj[start_dof+1];j++) {
1305       PetscInt dof = adjncy[j];
1306       if (!PetscBTLookup(touched,dof)) {
1307         ierr = PetscBTSet(touched,dof);CHKERRQ(ierr);
1308         queue_tip[n] = dof;
1309         n++;
1310       }
1311     }
1312   }
1313   *n_added = n;
1314   PetscFunctionReturn(0);
1315 }
1316