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