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