xref: /petsc/src/ksp/pc/impls/bddc/bddcschurs.c (revision 4c6709b31df94c6e6cd21f7a1b246810db8030e6)
134a97f8cSStefano Zampini #include <../src/ksp/pc/impls/bddc/bddc.h>
234a97f8cSStefano Zampini #include <../src/ksp/pc/impls/bddc/bddcprivate.h>
308122e43SStefano Zampini #include <petscblaslapack.h>
434a97f8cSStefano Zampini 
53202ece2SStefano Zampini PETSC_STATIC_INLINE PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt*,PetscInt,PetscBT,PetscInt*,PetscInt*,PetscInt*);
65ec10c6aSStefano Zampini static PetscErrorCode PCBDDCComputeExplicitSchur(Mat,PetscBool,MatReuse,Mat*);
73202ece2SStefano Zampini 
83202ece2SStefano Zampini #undef __FUNCT__
93202ece2SStefano Zampini #define __FUNCT__ "PCBDDCComputeExplicitSchur"
105ec10c6aSStefano Zampini static PetscErrorCode PCBDDCComputeExplicitSchur(Mat M, PetscBool issym, MatReuse reuse, Mat *S)
113202ece2SStefano Zampini {
123202ece2SStefano Zampini   Mat            B, C, D, Bd, Cd, AinvBd;
133202ece2SStefano Zampini   KSP            ksp;
143202ece2SStefano Zampini   PC             pc;
153202ece2SStefano Zampini   PetscBool      isLU, isILU, isCHOL, Bdense, Cdense;
163202ece2SStefano Zampini   PetscReal      fill = 2.0;
17f11841e3SStefano Zampini   PetscInt       n_I;
183202ece2SStefano Zampini   PetscMPIInt    size;
193202ece2SStefano Zampini   PetscErrorCode ierr;
203202ece2SStefano Zampini 
213202ece2SStefano Zampini   PetscFunctionBegin;
223202ece2SStefano Zampini   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)M),&size);CHKERRQ(ierr);
233202ece2SStefano Zampini   if (size != 1) {
243202ece2SStefano Zampini     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not for parallel matrices");
253202ece2SStefano Zampini   }
26f11841e3SStefano Zampini   if (reuse == MAT_REUSE_MATRIX) {
27f11841e3SStefano Zampini     PetscBool Sdense;
28f11841e3SStefano Zampini 
29f11841e3SStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)*S, MATSEQDENSE, &Sdense);CHKERRQ(ierr);
30f11841e3SStefano Zampini     if (!Sdense) {
31f11841e3SStefano Zampini       SETERRQ(PetscObjectComm((PetscObject)M),PETSC_ERR_SUP,"S should dense");
32f11841e3SStefano Zampini     }
33f11841e3SStefano Zampini   }
343202ece2SStefano Zampini   ierr = MatSchurComplementGetSubMatrices(M, NULL, NULL, &B, &C, &D);CHKERRQ(ierr);
353202ece2SStefano Zampini   ierr = MatSchurComplementGetKSP(M, &ksp);CHKERRQ(ierr);
363202ece2SStefano Zampini   ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr);
373202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) pc, PCLU, &isLU);CHKERRQ(ierr);
383202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) pc, PCILU, &isILU);CHKERRQ(ierr);
393202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) pc, PCCHOLESKY, &isCHOL);CHKERRQ(ierr);
403202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) B, MATSEQDENSE, &Bdense);CHKERRQ(ierr);
413202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) C, MATSEQDENSE, &Cdense);CHKERRQ(ierr);
42f11841e3SStefano Zampini   ierr = MatGetSize(B,&n_I,NULL);CHKERRQ(ierr);
43f11841e3SStefano Zampini   if (n_I) {
443202ece2SStefano Zampini     if (!Bdense) {
453202ece2SStefano Zampini       ierr = MatConvert(B, MATSEQDENSE, MAT_INITIAL_MATRIX, &Bd);CHKERRQ(ierr);
463202ece2SStefano Zampini     } else {
473202ece2SStefano Zampini       Bd = B;
483202ece2SStefano Zampini     }
493202ece2SStefano Zampini 
503202ece2SStefano Zampini     if (isLU || isILU || isCHOL) {
513202ece2SStefano Zampini       Mat fact;
523202ece2SStefano Zampini 
533202ece2SStefano Zampini       ierr = KSPSetUp(ksp);CHKERRQ(ierr);
543202ece2SStefano Zampini       ierr = PCFactorGetMatrix(pc, &fact);CHKERRQ(ierr);
553202ece2SStefano Zampini       ierr = MatDuplicate(Bd, MAT_DO_NOT_COPY_VALUES, &AinvBd);CHKERRQ(ierr);
563202ece2SStefano Zampini       ierr = MatMatSolve(fact, Bd, AinvBd);CHKERRQ(ierr);
573202ece2SStefano Zampini     } else {
583202ece2SStefano Zampini       Mat Ainvd;
593202ece2SStefano Zampini 
603202ece2SStefano Zampini       ierr = PCComputeExplicitOperator(pc, &Ainvd);CHKERRQ(ierr);
613202ece2SStefano Zampini       ierr = MatMatMult(Ainvd, Bd, MAT_INITIAL_MATRIX, fill, &AinvBd);CHKERRQ(ierr);
623202ece2SStefano Zampini       ierr = MatDestroy(&Ainvd);CHKERRQ(ierr);
633202ece2SStefano Zampini     }
645ec10c6aSStefano Zampini     if (!Bdense & !issym) {
653202ece2SStefano Zampini       ierr = MatDestroy(&Bd);CHKERRQ(ierr);
663202ece2SStefano Zampini     }
675ec10c6aSStefano Zampini 
685ec10c6aSStefano Zampini     if (!issym) {
693202ece2SStefano Zampini       if (!Cdense) {
703202ece2SStefano Zampini         ierr = MatConvert(C, MATSEQDENSE, MAT_INITIAL_MATRIX, &Cd);CHKERRQ(ierr);
713202ece2SStefano Zampini       } else {
723202ece2SStefano Zampini         Cd = C;
733202ece2SStefano Zampini       }
745ec10c6aSStefano Zampini       ierr = MatMatMult(Cd, AinvBd, reuse, fill, S);CHKERRQ(ierr);
753202ece2SStefano Zampini       if (!Cdense) {
763202ece2SStefano Zampini         ierr = MatDestroy(&Cd);CHKERRQ(ierr);
773202ece2SStefano Zampini       }
785ec10c6aSStefano Zampini     } else {
795ec10c6aSStefano Zampini       ierr = MatTransposeMatMult(Bd, AinvBd, reuse, fill, S);CHKERRQ(ierr);
805ec10c6aSStefano Zampini       if (!Bdense) {
815ec10c6aSStefano Zampini         ierr = MatDestroy(&Bd);CHKERRQ(ierr);
825ec10c6aSStefano Zampini       }
835ec10c6aSStefano Zampini     }
845ec10c6aSStefano Zampini     ierr = MatDestroy(&AinvBd);CHKERRQ(ierr);
85f11841e3SStefano Zampini   }
863202ece2SStefano Zampini 
873202ece2SStefano Zampini   if (D) {
883202ece2SStefano Zampini     Mat       Dd;
893202ece2SStefano Zampini     PetscBool Ddense;
903202ece2SStefano Zampini 
913202ece2SStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)D,MATSEQDENSE,&Ddense);CHKERRQ(ierr);
923202ece2SStefano Zampini     if (!Ddense) {
933202ece2SStefano Zampini       ierr = MatConvert(D, MATSEQDENSE, MAT_INITIAL_MATRIX, &Dd);CHKERRQ(ierr);
943202ece2SStefano Zampini     } else {
953202ece2SStefano Zampini       Dd = D;
963202ece2SStefano Zampini     }
97f11841e3SStefano Zampini     if (n_I) {
983202ece2SStefano Zampini       ierr = MatAYPX(*S,-1.0,Dd,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
99f11841e3SStefano Zampini     } else {
100f11841e3SStefano Zampini       if (reuse == MAT_INITIAL_MATRIX) {
101f11841e3SStefano Zampini         ierr = MatDuplicate(Dd,MAT_COPY_VALUES,S);CHKERRQ(ierr);
102f11841e3SStefano Zampini       } else {
103f11841e3SStefano Zampini         ierr = MatCopy(Dd,*S,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
104f11841e3SStefano Zampini       }
105f11841e3SStefano Zampini     }
1063202ece2SStefano Zampini     if (!Ddense) {
1073202ece2SStefano Zampini       ierr = MatDestroy(&Dd);CHKERRQ(ierr);
1083202ece2SStefano Zampini     }
1093202ece2SStefano Zampini   } else {
1103202ece2SStefano Zampini     ierr = MatScale(*S,-1.0);CHKERRQ(ierr);
1113202ece2SStefano Zampini   }
1123202ece2SStefano Zampini   PetscFunctionReturn(0);
1133202ece2SStefano Zampini }
11434a97f8cSStefano Zampini 
11534a97f8cSStefano Zampini #undef __FUNCT__
1161580ed26SStefano Zampini #define __FUNCT__ "PCBDDCSubSchursSetUp"
1179552c7c7SStefano Zampini PetscErrorCode PCBDDCSubSchursSetUp(PCBDDCSubSchurs sub_schurs, PetscInt xadj[], PetscInt adjncy[], PetscInt nlayers, PetscBool compute_Stilda, PetscBool deluxe, PetscBool invert_Stildas, PetscBool use_edges, PetscBool use_faces)
118b1b3d7a2SStefano Zampini {
119b1b3d7a2SStefano Zampini   Mat                    A_II,A_IB,A_BI,A_BB;
120b1b3d7a2SStefano Zampini   Mat                    AE_II,*AE_IE,*AE_EI,*AE_EE;
121883469d8SStefano Zampini   Mat                    S_all,global_schur_subsets,work_mat;
12208122e43SStefano Zampini   Mat                    S_Ej_tilda_all,S_Ej_inv_all;
1235db18549SStefano Zampini   ISLocalToGlobalMapping l2gmap_subsets;
1245db18549SStefano Zampini   IS                     is_I,*is_subset_B,temp_is;
125d648f858SStefano Zampini   PetscInt               *nnz,*all_local_idx_G,*all_local_idx_N;
1265ec10c6aSStefano Zampini   PetscInt               i,subset_size,max_subset_size_seq;
127883469d8SStefano Zampini   PetscInt               extra,local_size,global_size;
12808122e43SStefano Zampini   PetscBLASInt           B_N,B_ierr,B_lwork,*pivots;
12906a4b1faSStefano Zampini   PetscScalar            *Bwork;
130b1b3d7a2SStefano Zampini   PetscErrorCode         ierr;
131b1b3d7a2SStefano Zampini 
132b1b3d7a2SStefano Zampini   PetscFunctionBegin;
133b1b3d7a2SStefano Zampini   /* get Schur complement matrices */
134883469d8SStefano Zampini   if (!sub_schurs->use_mumps) {
135f11841e3SStefano Zampini     PetscBool isseqaij;
136f11841e3SStefano Zampini 
137b7eb3628SStefano Zampini     if (compute_Stilda) {
1383dc780c3SStefano Zampini       SETERRQ(PetscObjectComm((PetscObject)sub_schurs->l2gmap),PETSC_ERR_SUP,"Adaptive selection of constraints requires MUMPS");
139b7eb3628SStefano Zampini     }
140b1b3d7a2SStefano Zampini     ierr = MatSchurComplementGetSubMatrices(sub_schurs->S,&A_II,NULL,&A_IB,&A_BI,&A_BB);CHKERRQ(ierr);
141f11841e3SStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)A_BB,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
142f11841e3SStefano Zampini     if (!isseqaij) {
143f11841e3SStefano Zampini       ierr = MatConvert(A_BB,MATSEQAIJ,MAT_REUSE_MATRIX,&A_BB);CHKERRQ(ierr);
144f11841e3SStefano Zampini       ierr = MatConvert(A_IB,MATSEQAIJ,MAT_REUSE_MATRIX,&A_IB);CHKERRQ(ierr);
145f11841e3SStefano Zampini       ierr = MatConvert(A_BI,MATSEQAIJ,MAT_REUSE_MATRIX,&A_BI);CHKERRQ(ierr);
146f11841e3SStefano Zampini     }
147b1b3d7a2SStefano Zampini     ierr = PetscMalloc4(sub_schurs->n_subs,&is_subset_B,
148b1b3d7a2SStefano Zampini                         sub_schurs->n_subs,&AE_IE,
149b1b3d7a2SStefano Zampini                         sub_schurs->n_subs,&AE_EI,
150b1b3d7a2SStefano Zampini                         sub_schurs->n_subs,&AE_EE);CHKERRQ(ierr);
151a58a30b4SStefano Zampini   } else {
152a58a30b4SStefano Zampini     is_subset_B = NULL;
153a58a30b4SStefano Zampini     AE_IE = NULL;
154a58a30b4SStefano Zampini     AE_EI = NULL;
155a58a30b4SStefano Zampini     AE_EE = NULL;
156b1b3d7a2SStefano Zampini   }
157b1b3d7a2SStefano Zampini 
158b1b3d7a2SStefano Zampini   /* determine interior problems */
159b96c3477SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_I_layer);CHKERRQ(ierr);
1603dc780c3SStefano Zampini   ierr = ISGetLocalSize(sub_schurs->is_I,&i);CHKERRQ(ierr);
1613dc780c3SStefano Zampini   if (nlayers >= 0 && i) { /* Interior problems can be different from the original one */
162b1b3d7a2SStefano Zampini     PetscBT                touched;
163b1b3d7a2SStefano Zampini     const PetscInt*        idx_B;
164b1b3d7a2SStefano Zampini     PetscInt               n_I,n_B,n_local_dofs,n_prev_added,j,layer,*local_numbering;
165b1b3d7a2SStefano Zampini 
1663dc780c3SStefano Zampini     if (xadj == NULL || adjncy == NULL) {
1673dc780c3SStefano Zampini       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cannot request layering without adjacency");
1683dc780c3SStefano Zampini     }
169b1b3d7a2SStefano Zampini     /* get sizes */
170b1b3d7a2SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr);
171b1b3d7a2SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_B,&n_B);CHKERRQ(ierr);
172b1b3d7a2SStefano Zampini 
173b1b3d7a2SStefano Zampini     ierr = PetscMalloc1(n_I+n_B,&local_numbering);CHKERRQ(ierr);
174b1b3d7a2SStefano Zampini     ierr = PetscBTCreate(n_I+n_B,&touched);CHKERRQ(ierr);
175b1b3d7a2SStefano Zampini     ierr = PetscBTMemzero(n_I+n_B,touched);CHKERRQ(ierr);
176b1b3d7a2SStefano Zampini 
177b1b3d7a2SStefano Zampini     /* all boundary dofs must be skipped when adding layers */
178b1b3d7a2SStefano Zampini     ierr = ISGetIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr);
179b1b3d7a2SStefano Zampini     for (j=0;j<n_B;j++) {
180b1b3d7a2SStefano Zampini       ierr = PetscBTSet(touched,idx_B[j]);CHKERRQ(ierr);
181b1b3d7a2SStefano Zampini     }
182b1b3d7a2SStefano Zampini     ierr = PetscMemcpy(local_numbering,idx_B,n_B*sizeof(PetscInt));CHKERRQ(ierr);
183b1b3d7a2SStefano Zampini     ierr = ISRestoreIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr);
184b1b3d7a2SStefano Zampini 
185b1b3d7a2SStefano Zampini     /* add prescribed number of layers of dofs */
186b1b3d7a2SStefano Zampini     n_local_dofs = n_B;
187b1b3d7a2SStefano Zampini     n_prev_added = n_B;
188b1b3d7a2SStefano Zampini     for (layer=0;layer<nlayers;layer++) {
189b1b3d7a2SStefano Zampini       PetscInt n_added;
190b1b3d7a2SStefano Zampini       if (n_local_dofs == n_I+n_B) break;
191b1b3d7a2SStefano Zampini       if (n_local_dofs > n_I+n_B) {
192b1b3d7a2SStefano Zampini         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);
193b1b3d7a2SStefano Zampini       }
194b1b3d7a2SStefano Zampini       ierr = PCBDDCAdjGetNextLayer_Private(local_numbering+n_local_dofs,n_prev_added,touched,xadj,adjncy,&n_added);CHKERRQ(ierr);
195b1b3d7a2SStefano Zampini       n_prev_added = n_added;
196b1b3d7a2SStefano Zampini       n_local_dofs += n_added;
197b1b3d7a2SStefano Zampini       if (!n_added) break;
198b1b3d7a2SStefano Zampini     }
199b1b3d7a2SStefano Zampini     ierr = PetscBTDestroy(&touched);CHKERRQ(ierr);
200b1b3d7a2SStefano Zampini 
201883469d8SStefano Zampini     /* IS for I layer dofs in original numbering */
20268270318SStefano Zampini     ierr = ISCreateGeneral(PetscObjectComm((PetscObject)sub_schurs->is_I),n_local_dofs-n_B,local_numbering+n_B,PETSC_COPY_VALUES,&sub_schurs->is_I_layer);CHKERRQ(ierr);
203b1b3d7a2SStefano Zampini     ierr = PetscFree(local_numbering);CHKERRQ(ierr);
20468270318SStefano Zampini     ierr = ISSort(sub_schurs->is_I_layer);CHKERRQ(ierr);
205883469d8SStefano Zampini     /* IS for I layer dofs in I numbering */
206883469d8SStefano Zampini     if (!sub_schurs->use_mumps) {
207b1b3d7a2SStefano Zampini       ISLocalToGlobalMapping ItoNmap;
208b1b3d7a2SStefano Zampini       ierr = ISLocalToGlobalMappingCreateIS(sub_schurs->is_I,&ItoNmap);CHKERRQ(ierr);
20968270318SStefano Zampini       ierr = ISGlobalToLocalMappingApplyIS(ItoNmap,IS_GTOLM_DROP,sub_schurs->is_I_layer,&is_I);CHKERRQ(ierr);
210b1b3d7a2SStefano Zampini       ierr = ISLocalToGlobalMappingDestroy(&ItoNmap);CHKERRQ(ierr);
211b1b3d7a2SStefano Zampini 
212b1b3d7a2SStefano Zampini       /* II block */
213b1b3d7a2SStefano Zampini       ierr = MatGetSubMatrix(A_II,is_I,is_I,MAT_INITIAL_MATRIX,&AE_II);CHKERRQ(ierr);
214b1b3d7a2SStefano Zampini     }
215b1b3d7a2SStefano Zampini   } else {
216b1b3d7a2SStefano Zampini     PetscInt n_I;
217b1b3d7a2SStefano Zampini 
218b1b3d7a2SStefano Zampini     /* IS for I dofs in original numbering */
219b1b3d7a2SStefano Zampini     ierr = PetscObjectReference((PetscObject)sub_schurs->is_I);CHKERRQ(ierr);
22068270318SStefano Zampini     sub_schurs->is_I_layer = sub_schurs->is_I;
221b1b3d7a2SStefano Zampini 
222b1b3d7a2SStefano Zampini     /* IS for I dofs in I numbering (strided 1) */
223883469d8SStefano Zampini     if (!sub_schurs->use_mumps) {
224b1b3d7a2SStefano Zampini       ierr = ISGetSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr);
225b1b3d7a2SStefano Zampini       ierr = ISCreateStride(PetscObjectComm((PetscObject)sub_schurs->is_I),n_I,0,1,&is_I);CHKERRQ(ierr);
226b1b3d7a2SStefano Zampini 
227b1b3d7a2SStefano Zampini       /* II block is the same */
228b1b3d7a2SStefano Zampini       ierr = PetscObjectReference((PetscObject)A_II);CHKERRQ(ierr);
229b1b3d7a2SStefano Zampini       AE_II = A_II;
230b1b3d7a2SStefano Zampini     }
231b1b3d7a2SStefano Zampini   }
232883469d8SStefano Zampini   /* Get info on subset sizes and sum of all subsets sizes */
2335ec10c6aSStefano Zampini   max_subset_size_seq = 0;
234883469d8SStefano Zampini   local_size = 0;
235883469d8SStefano Zampini   for (i=0;i<sub_schurs->n_subs_seq;i++) {
236883469d8SStefano Zampini     PetscInt j = sub_schurs->index_sequential[i];
237b96c3477SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_subs[j],&subset_size);CHKERRQ(ierr);
2385ec10c6aSStefano Zampini     max_subset_size_seq = PetscMax(subset_size,max_subset_size_seq);
239883469d8SStefano Zampini     local_size += subset_size;
240883469d8SStefano Zampini   }
241883469d8SStefano Zampini 
242883469d8SStefano Zampini   /* Work arrays for local indices */
243883469d8SStefano Zampini   extra = 0;
244883469d8SStefano Zampini   if (sub_schurs->use_mumps) {
245883469d8SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_I_layer,&extra);CHKERRQ(ierr);
246883469d8SStefano Zampini   }
247883469d8SStefano Zampini   ierr = PetscMalloc1(local_size+extra,&all_local_idx_N);CHKERRQ(ierr);
248883469d8SStefano Zampini   if (extra) {
249883469d8SStefano Zampini     const PetscInt *idxs;
250883469d8SStefano Zampini     ierr = ISGetIndices(sub_schurs->is_I_layer,&idxs);CHKERRQ(ierr);
251883469d8SStefano Zampini     ierr = PetscMemcpy(all_local_idx_N,idxs,extra*sizeof(PetscInt));CHKERRQ(ierr);
252883469d8SStefano Zampini     ierr = ISRestoreIndices(sub_schurs->is_I_layer,&idxs);CHKERRQ(ierr);
253883469d8SStefano Zampini   }
254883469d8SStefano Zampini   ierr = PetscMalloc1(local_size,&nnz);CHKERRQ(ierr);
255883469d8SStefano Zampini 
256883469d8SStefano Zampini   /* Get local indices in local numbering */
257883469d8SStefano Zampini   local_size = 0;
258883469d8SStefano Zampini   for (i=0;i<sub_schurs->n_subs_seq;i++) {
259883469d8SStefano Zampini     PetscInt j;
260883469d8SStefano Zampini     const    PetscInt *idxs;
261883469d8SStefano Zampini 
262883469d8SStefano Zampini     PetscInt local_problem_index = sub_schurs->index_sequential[i];
263b96c3477SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_subs[local_problem_index],&subset_size);CHKERRQ(ierr);
264b96c3477SStefano Zampini     ierr = ISGetIndices(sub_schurs->is_subs[local_problem_index],&idxs);CHKERRQ(ierr);
265883469d8SStefano Zampini     /* subset indices in local numbering */
266883469d8SStefano Zampini     ierr = PetscMemcpy(all_local_idx_N+local_size+extra,idxs,subset_size*sizeof(PetscInt));CHKERRQ(ierr);
267b96c3477SStefano Zampini     ierr = ISRestoreIndices(sub_schurs->is_subs[local_problem_index],&idxs);CHKERRQ(ierr);
268883469d8SStefano Zampini     for (j=0;j<subset_size;j++) nnz[local_size+j] = subset_size;
269883469d8SStefano Zampini     local_size += subset_size;
270883469d8SStefano Zampini   }
271883469d8SStefano Zampini 
272883469d8SStefano Zampini   S_all = NULL;
2739ab7bb16SStefano Zampini   sub_schurs->is_hermitian = PETSC_FALSE;
2749ab7bb16SStefano Zampini   sub_schurs->is_posdef = PETSC_FALSE;
2759552c7c7SStefano Zampini   if (sub_schurs->A) {
2769ab7bb16SStefano Zampini     ierr = MatIsHermitian(sub_schurs->A,0.0,&sub_schurs->is_hermitian);CHKERRQ(ierr);
2779ab7bb16SStefano Zampini     if (sub_schurs->is_hermitian) {
2789ab7bb16SStefano Zampini       PetscScalar val;
2799ab7bb16SStefano Zampini       Vec         vec1,vec2;
2809ab7bb16SStefano Zampini 
2819ab7bb16SStefano Zampini       ierr = MatCreateVecs(sub_schurs->A,&vec1,&vec2);CHKERRQ(ierr);
2829ab7bb16SStefano Zampini       ierr = VecSetRandom(vec1,NULL);
2839ab7bb16SStefano Zampini       ierr = VecCopy(vec1,vec2);CHKERRQ(ierr);
2849ab7bb16SStefano Zampini       ierr = MatMult(sub_schurs->A,vec2,vec1);CHKERRQ(ierr);
2859ab7bb16SStefano Zampini       ierr = VecDot(vec1,vec2,&val);CHKERRQ(ierr);
2869ab7bb16SStefano Zampini       if (PetscRealPart(val) > 0. && PetscImaginaryPart(val) == 0.) sub_schurs->is_posdef = PETSC_TRUE;
2879ab7bb16SStefano Zampini       ierr = VecDestroy(&vec1);CHKERRQ(ierr);
2889ab7bb16SStefano Zampini       ierr = VecDestroy(&vec2);CHKERRQ(ierr);
2899ab7bb16SStefano Zampini     }
2909552c7c7SStefano Zampini   }
2912972d61bSStefano Zampini 
2923dc780c3SStefano Zampini   if (sub_schurs->n_subs) {
293883469d8SStefano Zampini     if (!sub_schurs->use_mumps) {
294883469d8SStefano Zampini       /* subsets in original and boundary numbering */
295883469d8SStefano Zampini       for (i=0;i<sub_schurs->n_subs;i++) {
296b96c3477SStefano Zampini         ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,sub_schurs->is_subs[i],&is_subset_B[i]);CHKERRQ(ierr);
297b1b3d7a2SStefano Zampini       }
298b1b3d7a2SStefano Zampini 
299b1b3d7a2SStefano Zampini       /* EE block */
300b1b3d7a2SStefano Zampini       for (i=0;i<sub_schurs->n_subs;i++) {
301b1b3d7a2SStefano Zampini         ierr = MatGetSubMatrix(A_BB,is_subset_B[i],is_subset_B[i],MAT_INITIAL_MATRIX,&AE_EE[i]);CHKERRQ(ierr);
302b1b3d7a2SStefano Zampini       }
303aa83b6aeSStefano Zampini 
304b1b3d7a2SStefano Zampini       /* IE block */
305b1b3d7a2SStefano Zampini       for (i=0;i<sub_schurs->n_subs;i++) {
306b1b3d7a2SStefano Zampini         ierr = MatGetSubMatrix(A_IB,is_I,is_subset_B[i],MAT_INITIAL_MATRIX,&AE_IE[i]);CHKERRQ(ierr);
307b1b3d7a2SStefano Zampini       }
308aa83b6aeSStefano Zampini 
309b1b3d7a2SStefano Zampini       /* EI block */
310b1b3d7a2SStefano Zampini       for (i=0;i<sub_schurs->n_subs;i++) {
311aa83b6aeSStefano Zampini         if (sub_schurs->is_hermitian) {
312aa83b6aeSStefano Zampini           ierr = MatCreateTranspose(AE_IE[i],&AE_EI[i]);CHKERRQ(ierr);
313aa83b6aeSStefano Zampini         } else {
314b1b3d7a2SStefano Zampini           ierr = MatGetSubMatrix(A_BI,is_subset_B[i],is_I,MAT_INITIAL_MATRIX,&AE_EI[i]);CHKERRQ(ierr);
315b1b3d7a2SStefano Zampini         }
316aa83b6aeSStefano Zampini       }
317b1b3d7a2SStefano Zampini 
318b1b3d7a2SStefano Zampini       /* setup Schur complements on subset */
319b1b3d7a2SStefano Zampini       for (i=0;i<sub_schurs->n_subs;i++) {
320b96c3477SStefano Zampini         ierr = MatDestroy(&sub_schurs->S_Ej[i]);CHKERRQ(ierr);
321b1b3d7a2SStefano Zampini         ierr = MatCreateSchurComplement(AE_II,AE_II,AE_IE[i],AE_EI[i],AE_EE[i],&sub_schurs->S_Ej[i]);CHKERRQ(ierr);
322a58a30b4SStefano Zampini         ierr = MatDestroy(&AE_EE[i]);CHKERRQ(ierr);
323a58a30b4SStefano Zampini         ierr = MatDestroy(&AE_IE[i]);CHKERRQ(ierr);
324a58a30b4SStefano Zampini         ierr = MatDestroy(&AE_EI[i]);CHKERRQ(ierr);
325b1b3d7a2SStefano Zampini         if (AE_II == A_II) { /* we can reuse the same ksp */
326b1b3d7a2SStefano Zampini           KSP ksp;
327b1b3d7a2SStefano Zampini           ierr = MatSchurComplementGetKSP(sub_schurs->S,&ksp);CHKERRQ(ierr);
328b1b3d7a2SStefano Zampini           ierr = MatSchurComplementSetKSP(sub_schurs->S_Ej[i],ksp);CHKERRQ(ierr);
329b1b3d7a2SStefano Zampini         } else { /* build new ksp object which inherits ksp and pc types from the original one */
330b1b3d7a2SStefano Zampini           KSP      origksp,schurksp;
331b1b3d7a2SStefano Zampini           PC       origpc,schurpc;
332b1b3d7a2SStefano Zampini           KSPType  ksp_type;
333b1b3d7a2SStefano Zampini           PCType   pc_type;
334b1b3d7a2SStefano Zampini           PetscInt n_internal;
335b1b3d7a2SStefano Zampini 
336b1b3d7a2SStefano Zampini           ierr = MatSchurComplementGetKSP(sub_schurs->S,&origksp);CHKERRQ(ierr);
337b1b3d7a2SStefano Zampini           ierr = MatSchurComplementGetKSP(sub_schurs->S_Ej[i],&schurksp);CHKERRQ(ierr);
338b1b3d7a2SStefano Zampini           ierr = KSPGetType(origksp,&ksp_type);CHKERRQ(ierr);
339b1b3d7a2SStefano Zampini           ierr = KSPSetType(schurksp,ksp_type);CHKERRQ(ierr);
340b1b3d7a2SStefano Zampini           ierr = KSPGetPC(schurksp,&schurpc);CHKERRQ(ierr);
341b1b3d7a2SStefano Zampini           ierr = KSPGetPC(origksp,&origpc);CHKERRQ(ierr);
342b1b3d7a2SStefano Zampini           ierr = PCGetType(origpc,&pc_type);CHKERRQ(ierr);
343b1b3d7a2SStefano Zampini           ierr = PCSetType(schurpc,pc_type);CHKERRQ(ierr);
344b1b3d7a2SStefano Zampini           ierr = ISGetSize(is_I,&n_internal);CHKERRQ(ierr);
345b1b3d7a2SStefano Zampini           if (n_internal) { /* UMFPACK gives error with 0 sized problems */
346b1b3d7a2SStefano Zampini             MatSolverPackage solver=NULL;
347b1b3d7a2SStefano Zampini             ierr = PCFactorGetMatSolverPackage(origpc,(const MatSolverPackage*)&solver);CHKERRQ(ierr);
348b1b3d7a2SStefano Zampini             if (solver) {
349b1b3d7a2SStefano Zampini               ierr = PCFactorSetMatSolverPackage(schurpc,solver);CHKERRQ(ierr);
350b1b3d7a2SStefano Zampini             }
351b1b3d7a2SStefano Zampini           }
352b1b3d7a2SStefano Zampini           ierr = KSPSetUp(schurksp);CHKERRQ(ierr);
353b1b3d7a2SStefano Zampini         }
354b1b3d7a2SStefano Zampini       }
355b1b3d7a2SStefano Zampini       /* free */
356b1b3d7a2SStefano Zampini       ierr = ISDestroy(&is_I);CHKERRQ(ierr);
357b1b3d7a2SStefano Zampini       ierr = MatDestroy(&AE_II);CHKERRQ(ierr);
358b1b3d7a2SStefano Zampini       for (i=0;i<sub_schurs->n_subs;i++) {
359b1b3d7a2SStefano Zampini         ierr = ISDestroy(&is_subset_B[i]);CHKERRQ(ierr);
360b1b3d7a2SStefano Zampini       }
361b1b3d7a2SStefano Zampini       ierr = PetscFree4(is_subset_B,AE_IE,AE_EI,AE_EE);CHKERRQ(ierr);
362883469d8SStefano Zampini #if defined(PETSC_HAVE_MUMPS)
363883469d8SStefano Zampini     } else {
364883469d8SStefano Zampini       Mat           A,F;
365883469d8SStefano Zampini       IS            is_A_all;
366883469d8SStefano Zampini       PetscInt      *idxs_schur,n_I;
367883469d8SStefano Zampini 
368883469d8SStefano Zampini       /* get working mat */
369883469d8SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_I_layer,&n_I);CHKERRQ(ierr);
370883469d8SStefano Zampini       ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size+n_I,all_local_idx_N,PETSC_COPY_VALUES,&is_A_all);CHKERRQ(ierr);
371d648f858SStefano Zampini       ierr = MatGetSubMatrixUnsorted(sub_schurs->A,is_A_all,is_A_all,&A);CHKERRQ(ierr);
372883469d8SStefano Zampini       ierr = ISDestroy(&is_A_all);CHKERRQ(ierr);
373883469d8SStefano Zampini 
37408122e43SStefano Zampini       if (n_I) {
3759ab7bb16SStefano Zampini         if (sub_schurs->is_hermitian && sub_schurs->is_posdef) {
376883469d8SStefano Zampini           ierr = MatGetFactor(A,MATSOLVERMUMPS,MAT_FACTOR_CHOLESKY,&F);CHKERRQ(ierr);
377883469d8SStefano Zampini         } else {
378883469d8SStefano Zampini           ierr = MatGetFactor(A,MATSOLVERMUMPS,MAT_FACTOR_LU,&F);CHKERRQ(ierr);
379883469d8SStefano Zampini         }
380883469d8SStefano Zampini 
381883469d8SStefano Zampini         /* subsets ordered last */
382883469d8SStefano Zampini         ierr = PetscMalloc1(local_size,&idxs_schur);CHKERRQ(ierr);
383883469d8SStefano Zampini         for (i=0;i<local_size;i++) {
384883469d8SStefano Zampini           idxs_schur[i] = n_I+i+1;
385883469d8SStefano Zampini         }
386883469d8SStefano Zampini         ierr = MatMumpsSetSchurIndices(F,local_size,idxs_schur);CHKERRQ(ierr);
387883469d8SStefano Zampini         ierr = PetscFree(idxs_schur);CHKERRQ(ierr);
388883469d8SStefano Zampini 
389883469d8SStefano Zampini         /* factorization step */
3909ab7bb16SStefano Zampini         if (sub_schurs->is_hermitian && sub_schurs->is_posdef) {
391883469d8SStefano Zampini           ierr = MatCholeskyFactorSymbolic(F,A,NULL,NULL);CHKERRQ(ierr);
392883469d8SStefano Zampini           ierr = MatCholeskyFactorNumeric(F,A,NULL);CHKERRQ(ierr);
393883469d8SStefano Zampini         } else {
394883469d8SStefano Zampini           ierr = MatLUFactorSymbolic(F,A,NULL,NULL,NULL);CHKERRQ(ierr);
395883469d8SStefano Zampini           ierr = MatLUFactorNumeric(F,A,NULL);CHKERRQ(ierr);
396883469d8SStefano Zampini         }
397883469d8SStefano Zampini 
398883469d8SStefano Zampini         /* get explicit Schur Complement computed during numeric factorization */
399883469d8SStefano Zampini         ierr = MatMumpsGetSchurComplement(F,&S_all);CHKERRQ(ierr);
400883469d8SStefano Zampini         ierr = MatDestroy(&F);CHKERRQ(ierr);
40108122e43SStefano Zampini       } else {
40208122e43SStefano Zampini         ierr = MatConvert(A,MATSEQDENSE,MAT_INITIAL_MATRIX,&S_all);CHKERRQ(ierr);
40308122e43SStefano Zampini       }
40408122e43SStefano Zampini       ierr = MatDestroy(&A);CHKERRQ(ierr);
405883469d8SStefano Zampini #endif
406b1b3d7a2SStefano Zampini     }
4073dc780c3SStefano Zampini   } else {
4083dc780c3SStefano Zampini     ierr = PetscFree(nnz);CHKERRQ(ierr);
4093dc780c3SStefano Zampini     ierr = PetscFree(all_local_idx_N);CHKERRQ(ierr);
4103dc780c3SStefano Zampini   }
4115db18549SStefano Zampini 
4125db18549SStefano Zampini   if (!sub_schurs->n_subs_seq_g) {
4135db18549SStefano Zampini     PetscFunctionReturn(0);
4145db18549SStefano Zampini   }
4155db18549SStefano Zampini 
4165db18549SStefano Zampini   /* subset indices in local boundary numbering */
417b96c3477SStefano Zampini   if (!sub_schurs->is_Ej_all) {
418d648f858SStefano Zampini     PetscInt *all_local_idx_B;
419d648f858SStefano Zampini 
420d648f858SStefano Zampini     ierr = PetscMalloc1(local_size,&all_local_idx_B);CHKERRQ(ierr);
421883469d8SStefano Zampini     ierr = ISGlobalToLocalMappingApply(sub_schurs->BtoNmap,IS_GTOLM_DROP,local_size,all_local_idx_N+extra,&subset_size,all_local_idx_B);CHKERRQ(ierr);
4225db18549SStefano Zampini     if (subset_size != local_size) {
4235db18549SStefano Zampini       SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in sub_schurs serial (BtoNmap)! %d != %d\n",subset_size,local_size);
4245db18549SStefano Zampini     }
4255db18549SStefano Zampini     ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size,all_local_idx_B,PETSC_OWN_POINTER,&sub_schurs->is_Ej_all);CHKERRQ(ierr);
426b96c3477SStefano Zampini   }
4275db18549SStefano Zampini 
4283202ece2SStefano Zampini   /* Local matrix of all local Schur on subsets transposed */
429b96c3477SStefano Zampini   if (!sub_schurs->S_Ej_all) {
4305db18549SStefano Zampini     ierr = MatCreate(PETSC_COMM_SELF,&sub_schurs->S_Ej_all);CHKERRQ(ierr);
4315db18549SStefano Zampini     ierr = MatSetSizes(sub_schurs->S_Ej_all,PETSC_DECIDE,PETSC_DECIDE,local_size,local_size);CHKERRQ(ierr);
4325db18549SStefano Zampini     ierr = MatSetType(sub_schurs->S_Ej_all,MATAIJ);CHKERRQ(ierr);
4335db18549SStefano Zampini     ierr = MatSeqAIJSetPreallocation(sub_schurs->S_Ej_all,0,nnz);CHKERRQ(ierr);
434b96c3477SStefano Zampini   } else {
435b96c3477SStefano Zampini     ierr = MatZeroEntries(sub_schurs->S_Ej_all);CHKERRQ(ierr);
436b96c3477SStefano Zampini   }
437a1337663SStefano Zampini 
438a1337663SStefano Zampini   S_Ej_tilda_all = 0;
43908122e43SStefano Zampini   S_Ej_inv_all = 0;
44006a4b1faSStefano Zampini   Bwork = NULL;
441a58a30b4SStefano Zampini   pivots = NULL;
4422972d61bSStefano Zampini   if (sub_schurs->n_subs && deluxe && compute_Stilda && !sub_schurs->is_hermitian) { /* workspace needed only for GETRI */
44308122e43SStefano Zampini     PetscScalar lwork;
44408122e43SStefano Zampini 
44508122e43SStefano Zampini     B_lwork = -1;
4465ec10c6aSStefano Zampini     ierr = PetscBLASIntCast(max_subset_size_seq,&B_N);CHKERRQ(ierr);
44708122e43SStefano Zampini     ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
44806a4b1faSStefano Zampini     PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,Bwork,&B_N,pivots,&lwork,&B_lwork,&B_ierr));
44908122e43SStefano Zampini     ierr = PetscFPTrapPop();CHKERRQ(ierr);
45008122e43SStefano Zampini     if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to GETRI Lapack routine %d",(int)B_ierr);
451*4c6709b3SStefano Zampini     ierr = PetscBLASIntCast((PetscInt)PetscRealPart(lwork),&B_lwork);CHKERRQ(ierr);
45206a4b1faSStefano Zampini     ierr = PetscMalloc2(B_lwork,&Bwork,max_subset_size_seq,&pivots);CHKERRQ(ierr);
45308122e43SStefano Zampini   }
45408122e43SStefano Zampini 
45508122e43SStefano Zampini   ierr = PetscBTMemzero(sub_schurs->n_subs,sub_schurs->computed_Stilda_subs);CHKERRQ(ierr);
456883469d8SStefano Zampini   if (!sub_schurs->use_mumps) {
4575ec10c6aSStefano Zampini     if (sub_schurs->n_subs_seq) {
4585ec10c6aSStefano Zampini       PetscScalar *work;
4595db18549SStefano Zampini       PetscInt    *dummy_idx;
4605db18549SStefano Zampini 
4615db18549SStefano Zampini       /* Work arrays */
4625ec10c6aSStefano Zampini       ierr = PetscMalloc2(max_subset_size_seq,&dummy_idx,max_subset_size_seq*max_subset_size_seq,&work);CHKERRQ(ierr);
4635ec10c6aSStefano Zampini 
4643202ece2SStefano Zampini       /* Loop on local problems end compute Schur complements explicitly */
4655db18549SStefano Zampini       local_size = 0;
4665db18549SStefano Zampini       for (i=0;i<sub_schurs->n_subs_seq;i++) {
4673202ece2SStefano Zampini         Mat       S_Ej_expl;
4683202ece2SStefano Zampini         PetscInt  j,lpi = sub_schurs->index_sequential[i];
4693202ece2SStefano Zampini         PetscBool Sdense;
4705db18549SStefano Zampini 
471b96c3477SStefano Zampini         ierr = ISGetLocalSize(sub_schurs->is_subs[lpi],&subset_size);CHKERRQ(ierr);
4725ec10c6aSStefano Zampini         ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&S_Ej_expl);CHKERRQ(ierr);
473aa83b6aeSStefano Zampini         ierr = PCBDDCComputeExplicitSchur(sub_schurs->S_Ej[lpi],sub_schurs->is_hermitian,MAT_REUSE_MATRIX,&S_Ej_expl);CHKERRQ(ierr);
4743202ece2SStefano Zampini         ierr = PetscObjectTypeCompare((PetscObject)S_Ej_expl,MATSEQDENSE,&Sdense);CHKERRQ(ierr);
4753202ece2SStefano Zampini         if (Sdense) {
4765db18549SStefano Zampini           for (j=0;j<subset_size;j++) {
4775db18549SStefano Zampini             dummy_idx[j]=local_size+j;
4785db18549SStefano Zampini           }
4795ec10c6aSStefano Zampini           ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
4803202ece2SStefano Zampini         } else {
4813202ece2SStefano Zampini           SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not yet implemented for sparse matrices");
4825db18549SStefano Zampini         }
4833202ece2SStefano Zampini         local_size += subset_size;
48465d8bf0aSStefano Zampini         ierr = MatDestroy(&sub_schurs->S_Ej[lpi]);CHKERRQ(ierr);
48565d8bf0aSStefano Zampini         sub_schurs->S_Ej[lpi] = S_Ej_expl;
4863202ece2SStefano Zampini       }
4875ec10c6aSStefano Zampini       ierr = PetscFree2(dummy_idx,work);CHKERRQ(ierr);
4885ec10c6aSStefano Zampini     }
489883469d8SStefano Zampini   } else {
490883469d8SStefano Zampini     PetscInt    *dummy_idx,n_all;
4915ec10c6aSStefano Zampini     PetscScalar *work;
492883469d8SStefano Zampini 
493a1337663SStefano Zampini     if (compute_Stilda) {
494a1337663SStefano Zampini       ierr = MatCreate(PETSC_COMM_SELF,&S_Ej_tilda_all);CHKERRQ(ierr);
495a1337663SStefano Zampini       ierr = MatSetSizes(S_Ej_tilda_all,PETSC_DECIDE,PETSC_DECIDE,local_size,local_size);CHKERRQ(ierr);
496a1337663SStefano Zampini       ierr = MatSetType(S_Ej_tilda_all,MATAIJ);CHKERRQ(ierr);
497a1337663SStefano Zampini       ierr = MatSeqAIJSetPreallocation(S_Ej_tilda_all,0,nnz);CHKERRQ(ierr);
49808122e43SStefano Zampini       if (deluxe) {
49908122e43SStefano Zampini         ierr = MatCreate(PETSC_COMM_SELF,&S_Ej_inv_all);CHKERRQ(ierr);
50008122e43SStefano Zampini         ierr = MatSetSizes(S_Ej_inv_all,PETSC_DECIDE,PETSC_DECIDE,local_size,local_size);CHKERRQ(ierr);
50108122e43SStefano Zampini         ierr = MatSetType(S_Ej_inv_all,MATAIJ);CHKERRQ(ierr);
50208122e43SStefano Zampini         ierr = MatSeqAIJSetPreallocation(S_Ej_inv_all,0,nnz);CHKERRQ(ierr);
503a1337663SStefano Zampini       }
50408122e43SStefano Zampini     }
50508122e43SStefano Zampini 
5068275c6a8SStefano Zampini     if (S_all) { /* multilevel */
50708122e43SStefano Zampini       ierr = MatGetSize(S_all,&n_all,NULL);CHKERRQ(ierr);
5088275c6a8SStefano Zampini     } else {
5098275c6a8SStefano Zampini       n_all = 0;
5108275c6a8SStefano Zampini     }
51108122e43SStefano Zampini     local_size = 0;
51206a4b1faSStefano Zampini 
513883469d8SStefano Zampini     /* Work arrays */
51406a4b1faSStefano Zampini     if (compute_Stilda) {
51506a4b1faSStefano Zampini       ierr = PetscMalloc2(max_subset_size_seq,&dummy_idx,n_all*n_all,&work);CHKERRQ(ierr);
51606a4b1faSStefano Zampini     } else {
51706a4b1faSStefano Zampini       ierr = PetscMalloc2(max_subset_size_seq,&dummy_idx,0,&work);CHKERRQ(ierr);
51806a4b1faSStefano Zampini     }
519883469d8SStefano Zampini 
52065d8bf0aSStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
52165d8bf0aSStefano Zampini       IS  is_E;
52265d8bf0aSStefano Zampini       PetscScalar *vals;
52365d8bf0aSStefano Zampini       PetscInt j;
52465d8bf0aSStefano Zampini 
525b96c3477SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
526a1337663SStefano Zampini       for (j=0;j<subset_size;j++) {
527a1337663SStefano Zampini         dummy_idx[j]=local_size+j;
528a1337663SStefano Zampini       }
52965d8bf0aSStefano Zampini       ierr = ISCreateStride(PetscObjectComm((PetscObject)sub_schurs->is_I),subset_size,local_size,1,&is_E);CHKERRQ(ierr);
530b96c3477SStefano Zampini       if (sub_schurs->S_Ej[i]) {
531b96c3477SStefano Zampini         ierr = MatGetSubMatrix(S_all,is_E,is_E,MAT_REUSE_MATRIX,&sub_schurs->S_Ej[i]);CHKERRQ(ierr);
532b96c3477SStefano Zampini       } else {
533a1337663SStefano Zampini         ierr = MatGetSubMatrix(S_all,is_E,is_E,MAT_INITIAL_MATRIX,&sub_schurs->S_Ej[i]);CHKERRQ(ierr);
534b96c3477SStefano Zampini       }
535a1337663SStefano Zampini       ierr = MatDenseGetArray(sub_schurs->S_Ej[i],&vals);CHKERRQ(ierr);
536a1337663SStefano Zampini       ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,vals,INSERT_VALUES);CHKERRQ(ierr);
537a1337663SStefano Zampini       ierr = MatDenseRestoreArray(sub_schurs->S_Ej[i],&vals);CHKERRQ(ierr);
538a1337663SStefano Zampini 
5396da91a30SStefano Zampini       if (compute_Stilda && ((PetscBTLookup(sub_schurs->is_edge,i) && use_edges) || use_faces)) {
54065d8bf0aSStefano Zampini         Mat Stilda;
54106a4b1faSStefano Zampini 
54206a4b1faSStefano Zampini         ierr = PetscBTSet(sub_schurs->computed_Stilda_subs,i);CHKERRQ(ierr);
54365d8bf0aSStefano Zampini         if (sub_schurs->n_subs == 1) {
54465d8bf0aSStefano Zampini           ierr = PetscObjectReference((PetscObject)sub_schurs->S_Ej[i]);CHKERRQ(ierr);
54565d8bf0aSStefano Zampini           Stilda = sub_schurs->S_Ej[i];
54665d8bf0aSStefano Zampini         } else {
54765d8bf0aSStefano Zampini           KSP         ksp;
54865d8bf0aSStefano Zampini           PC          pc;
54965d8bf0aSStefano Zampini           Mat         S_EF,S_FE,S_FF,Stilda_impl;
55065d8bf0aSStefano Zampini           IS          is_F;
55106a4b1faSStefano Zampini           PetscInt    nc,cum;
5522b973097SStefano Zampini           PetscBool   chop=PETSC_FALSE;
55365d8bf0aSStefano Zampini 
55465d8bf0aSStefano Zampini           ierr = ISComplement(is_E,0,n_all,&is_F);CHKERRQ(ierr);
55506a4b1faSStefano Zampini           nc   = n_all - subset_size;
55606a4b1faSStefano Zampini           cum  = 0;
55706a4b1faSStefano Zampini           ierr = MatCreateSeqDense(PETSC_COMM_SELF,nc,nc,work+cum,&S_FF);CHKERRQ(ierr);
5582972d61bSStefano Zampini           cum += nc*nc;
55906a4b1faSStefano Zampini           ierr = MatCreateSeqDense(PETSC_COMM_SELF,nc,subset_size,work+cum,&S_FE);CHKERRQ(ierr);
5602972d61bSStefano Zampini           cum += nc*subset_size;
56106a4b1faSStefano Zampini           ierr = MatGetSubMatrix(S_all,is_F,is_F,MAT_REUSE_MATRIX,&S_FF);CHKERRQ(ierr);
56206a4b1faSStefano Zampini           ierr = MatGetSubMatrix(S_all,is_F,is_E,MAT_REUSE_MATRIX,&S_FE);CHKERRQ(ierr);
5635ec10c6aSStefano Zampini           if (sub_schurs->is_hermitian) { /* just need a placeholder for S_EF */
56406a4b1faSStefano Zampini             ierr = MatCreateTranspose(S_FE,&S_EF);CHKERRQ(ierr);
5655ec10c6aSStefano Zampini           } else {
5662972d61bSStefano Zampini             ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,nc,work+cum,&S_EF);CHKERRQ(ierr);
5672972d61bSStefano Zampini             cum += nc*subset_size;
56806a4b1faSStefano Zampini             ierr = MatGetSubMatrix(S_all,is_E,is_F,MAT_REUSE_MATRIX,&S_EF);CHKERRQ(ierr);
5695ec10c6aSStefano Zampini           }
570a1337663SStefano Zampini           ierr = ISDestroy(&is_F);CHKERRQ(ierr);
5712b973097SStefano Zampini           if (chop) {
572*4c6709b3SStefano Zampini             PetscReal eps=1.e-8;
5732b973097SStefano Zampini             ierr = MatChop(S_FF,eps);CHKERRQ(ierr);
5742b973097SStefano Zampini             ierr = MatConvert(S_FF,MATAIJ,MAT_REUSE_MATRIX,&S_FF);CHKERRQ(ierr);
5752b973097SStefano Zampini           }
576a1337663SStefano Zampini           ierr = MatCreateSchurComplement(S_FF,S_FF,S_FE,S_EF,sub_schurs->S_Ej[i],&Stilda_impl);CHKERRQ(ierr);
577a1337663SStefano Zampini           ierr = MatDestroy(&S_FF);CHKERRQ(ierr);
578a1337663SStefano Zampini           ierr = MatDestroy(&S_FE);CHKERRQ(ierr);
579a1337663SStefano Zampini           ierr = MatDestroy(&S_EF);CHKERRQ(ierr);
58065d8bf0aSStefano Zampini           ierr = MatSchurComplementGetKSP(Stilda_impl,&ksp);CHKERRQ(ierr);
58165d8bf0aSStefano Zampini           ierr = KSPSetType(ksp,KSPPREONLY);CHKERRQ(ierr);
58265d8bf0aSStefano Zampini           ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr);
5839ab7bb16SStefano Zampini           if (sub_schurs->is_hermitian && sub_schurs->is_posdef) {
58465d8bf0aSStefano Zampini             ierr = PCSetType(pc,PCCHOLESKY);CHKERRQ(ierr);
58565d8bf0aSStefano Zampini           } else {
58665d8bf0aSStefano Zampini             ierr = PCSetType(pc,PCLU);CHKERRQ(ierr);
58765d8bf0aSStefano Zampini           }
5882b973097SStefano Zampini           if (!chop) {
58965d8bf0aSStefano Zampini             ierr = PCFactorSetUseInPlace(pc,PETSC_TRUE);CHKERRQ(ierr);
5902b973097SStefano Zampini           } else {
5912b973097SStefano Zampini             ierr = PCFactorSetMatSolverPackage(pc,MATSOLVERMUMPS);CHKERRQ(ierr);
5922b973097SStefano Zampini           }
59365d8bf0aSStefano Zampini           ierr = KSPSetUp(ksp);CHKERRQ(ierr);
59406a4b1faSStefano Zampini           ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work+cum,&Stilda);CHKERRQ(ierr);
5955ec10c6aSStefano Zampini           ierr = PCBDDCComputeExplicitSchur(Stilda_impl,sub_schurs->is_hermitian,MAT_REUSE_MATRIX,&Stilda);CHKERRQ(ierr);
59665d8bf0aSStefano Zampini           ierr = MatDestroy(&Stilda_impl);CHKERRQ(ierr);
59765d8bf0aSStefano Zampini         }
5982b973097SStefano Zampini /*
5992b973097SStefano Zampini         PetscViewerPushFormat(PETSC_VIEWER_STDOUT_SELF,PETSC_VIEWER_ASCII_MATLAB);
6002b973097SStefano Zampini         ierr = MatView(Stilda,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
6012b973097SStefano Zampini         PetscViewerPopFormat(PETSC_VIEWER_STDOUT_SELF);
6022b973097SStefano Zampini */
603b76a6013SStefano Zampini         ierr = MatDenseGetArray(Stilda,&vals);CHKERRQ(ierr);
60408122e43SStefano Zampini         if (deluxe) { /* when using deluxe scaling, we need (S_1^-1+S_2^-1)^-1 */
60508122e43SStefano Zampini           PetscScalar *vals2;
60608122e43SStefano Zampini 
60708122e43SStefano Zampini           ierr = MatDenseGetArray(sub_schurs->S_Ej[i],&vals2);CHKERRQ(ierr);
60808122e43SStefano Zampini           ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr);
60908122e43SStefano Zampini           ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
6109552c7c7SStefano Zampini           if (invert_Stildas) { /* Stildas can be singular */
6112972d61bSStefano Zampini             if (!sub_schurs->is_hermitian) {
61208122e43SStefano Zampini               PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,vals,&B_N,pivots,&B_ierr));
61308122e43SStefano Zampini               if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
61406a4b1faSStefano Zampini               PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,vals,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
61508122e43SStefano Zampini               if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
6162972d61bSStefano Zampini             } else {
6172972d61bSStefano Zampini               PetscInt j,k;
6182972d61bSStefano Zampini 
6192972d61bSStefano Zampini               PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,vals,&B_N,&B_ierr));
6202972d61bSStefano Zampini               if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
6212972d61bSStefano Zampini               PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,vals,&B_N,&B_ierr));
6222972d61bSStefano Zampini               if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
6232972d61bSStefano Zampini               for (j=0;j<B_N;j++) {
6242972d61bSStefano Zampini                 for (k=j+1;k<B_N;k++) {
6252972d61bSStefano Zampini                   vals[k*B_N+j] = vals[j*B_N+k];
6269552c7c7SStefano Zampini                 }
6272972d61bSStefano Zampini               }
6282972d61bSStefano Zampini             }
6292972d61bSStefano Zampini           }
6302972d61bSStefano Zampini           if (!sub_schurs->is_hermitian) {
63108122e43SStefano Zampini             PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,vals2,&B_N,pivots,&B_ierr));
63208122e43SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
63306a4b1faSStefano Zampini             PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,vals2,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
63408122e43SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
6352972d61bSStefano Zampini           } else {
6362972d61bSStefano Zampini             PetscInt j,k;
6372972d61bSStefano Zampini 
6382972d61bSStefano Zampini             PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,vals2,&B_N,&B_ierr));
6392972d61bSStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
6402972d61bSStefano Zampini             PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,vals2,&B_N,&B_ierr));
6412972d61bSStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
6422972d61bSStefano Zampini             for (j=0;j<B_N;j++) {
6432972d61bSStefano Zampini               for (k=j+1;k<B_N;k++) {
6442972d61bSStefano Zampini                 vals2[k*B_N+j] = vals2[j*B_N+k];
6452972d61bSStefano Zampini               }
6462972d61bSStefano Zampini             }
6472972d61bSStefano Zampini           }
64808122e43SStefano Zampini           ierr = PetscFPTrapPop();CHKERRQ(ierr);
64908122e43SStefano Zampini           ierr = MatSetValues(S_Ej_inv_all,subset_size,dummy_idx,subset_size,dummy_idx,vals2,INSERT_VALUES);CHKERRQ(ierr);
65008122e43SStefano Zampini           ierr = MatDenseRestoreArray(sub_schurs->S_Ej[i],&vals2);CHKERRQ(ierr);
65108122e43SStefano Zampini         }
652a1337663SStefano Zampini         ierr = MatSetValues(S_Ej_tilda_all,subset_size,dummy_idx,subset_size,dummy_idx,vals,INSERT_VALUES);CHKERRQ(ierr);
653b76a6013SStefano Zampini         ierr = MatDenseRestoreArray(Stilda,&vals);CHKERRQ(ierr);
654b76a6013SStefano Zampini         ierr = MatDestroy(&Stilda);CHKERRQ(ierr);
65565d8bf0aSStefano Zampini       }
65665d8bf0aSStefano Zampini       ierr = ISDestroy(&is_E);CHKERRQ(ierr);
657883469d8SStefano Zampini       local_size += subset_size;
658883469d8SStefano Zampini     }
6595ec10c6aSStefano Zampini     ierr = PetscFree2(dummy_idx,work);CHKERRQ(ierr);
6605db18549SStefano Zampini   }
661a1337663SStefano Zampini   ierr = PetscFree(nnz);CHKERRQ(ierr);
662a1337663SStefano Zampini   ierr = MatDestroy(&S_all);CHKERRQ(ierr);
6635db18549SStefano Zampini   ierr = MatAssemblyBegin(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
6645db18549SStefano Zampini   ierr = MatAssemblyEnd(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
6658275c6a8SStefano Zampini   if (S_Ej_tilda_all) {
666a1337663SStefano Zampini     ierr = MatAssemblyBegin(S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
667a1337663SStefano Zampini     ierr = MatAssemblyEnd(S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
6688275c6a8SStefano Zampini   }
6698275c6a8SStefano Zampini   if (S_Ej_inv_all) {
67008122e43SStefano Zampini     ierr = MatAssemblyBegin(S_Ej_inv_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
67108122e43SStefano Zampini     ierr = MatAssemblyEnd(S_Ej_inv_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
67208122e43SStefano Zampini   }
673a1337663SStefano Zampini 
6745db18549SStefano Zampini   /* Global matrix of all assembled Schur on subsets */
675883469d8SStefano Zampini   ierr = PCBDDCSubsetNumbering(PetscObjectComm((PetscObject)sub_schurs->l2gmap),sub_schurs->l2gmap,local_size,all_local_idx_N+extra,PETSC_NULL,&global_size,&all_local_idx_G);CHKERRQ(ierr);
6765db18549SStefano Zampini   ierr = ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)sub_schurs->l2gmap),1,local_size,all_local_idx_G,PETSC_COPY_VALUES,&l2gmap_subsets);CHKERRQ(ierr);
6775db18549SStefano Zampini   ierr = MatCreateIS(PetscObjectComm((PetscObject)sub_schurs->l2gmap),1,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size,l2gmap_subsets,&work_mat);CHKERRQ(ierr);
6785db18549SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&l2gmap_subsets);CHKERRQ(ierr);
6795db18549SStefano Zampini   ierr = MatISSetLocalMat(work_mat,sub_schurs->S_Ej_all);CHKERRQ(ierr);
6803927de2eSStefano Zampini   ierr = MatCreate(PetscObjectComm((PetscObject)work_mat),&global_schur_subsets);CHKERRQ(ierr);
6813927de2eSStefano Zampini   ierr = MatSetSizes(global_schur_subsets,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size);CHKERRQ(ierr);
6823927de2eSStefano Zampini   ierr = MatSetType(global_schur_subsets,MATMPIAIJ);CHKERRQ(ierr);
6833927de2eSStefano Zampini   ierr = MatISSetMPIXAIJPreallocation_Private(work_mat,global_schur_subsets,PETSC_TRUE);CHKERRQ(ierr);
6843927de2eSStefano Zampini   ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
6855db18549SStefano Zampini   /* Get local part of (\sum_j S_Ej) */
686883469d8SStefano Zampini   ierr = PetscFree(all_local_idx_N);CHKERRQ(ierr);
687883469d8SStefano Zampini   ierr = ISCreateGeneral(PetscObjectComm((PetscObject)sub_schurs->l2gmap),local_size,all_local_idx_G,PETSC_OWN_POINTER,&temp_is);CHKERRQ(ierr);
688d648f858SStefano Zampini   ierr = MatDestroy(&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
689d648f858SStefano Zampini   ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
69008122e43SStefano Zampini 
6918275c6a8SStefano Zampini   if (S_Ej_tilda_all) {
692a1337663SStefano Zampini     ierr = MatISSetLocalMat(work_mat,S_Ej_tilda_all);CHKERRQ(ierr);
693a1337663SStefano Zampini     ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
694d648f858SStefano Zampini     ierr = MatDestroy(&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
695d648f858SStefano Zampini     ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
6968275c6a8SStefano Zampini   }
6978275c6a8SStefano Zampini   if (S_Ej_inv_all) {
69808122e43SStefano Zampini     PetscInt    cum;
69908122e43SStefano Zampini     PetscScalar *array,*array2;
70008122e43SStefano Zampini     ierr = MatISSetLocalMat(work_mat,S_Ej_inv_all);CHKERRQ(ierr);
70108122e43SStefano Zampini     ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
702d648f858SStefano Zampini     ierr = MatDestroy(&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
703d648f858SStefano Zampini     ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
70408122e43SStefano Zampini     /* invert blocks of sum_S_Ej_inv_all */
70508122e43SStefano Zampini     ierr = MatSeqAIJGetArray(sub_schurs->sum_S_Ej_inv_all,&array);CHKERRQ(ierr);
70608122e43SStefano Zampini     ierr = MatSeqAIJGetArray(sub_schurs->sum_S_Ej_tilda_all,&array2);CHKERRQ(ierr);
70708122e43SStefano Zampini     cum = 0;
70808122e43SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
70908122e43SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
71008122e43SStefano Zampini       if (PetscBTLookup(sub_schurs->computed_Stilda_subs,i)) {
71108122e43SStefano Zampini         ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr);
71208122e43SStefano Zampini         ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
7132972d61bSStefano Zampini         if (!sub_schurs->is_hermitian) {
71408122e43SStefano Zampini           PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,array+cum,&B_N,pivots,&B_ierr));
71508122e43SStefano Zampini           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
71606a4b1faSStefano Zampini           PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,array+cum,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
71708122e43SStefano Zampini           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
7182972d61bSStefano Zampini         } else {
7192972d61bSStefano Zampini           PetscInt j,k;
7202972d61bSStefano Zampini 
7212972d61bSStefano Zampini           PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,array+cum,&B_N,&B_ierr));
7222972d61bSStefano Zampini           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
7232972d61bSStefano Zampini           PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,array+cum,&B_N,&B_ierr));
7242972d61bSStefano Zampini           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
7252972d61bSStefano Zampini           for (j=0;j<B_N;j++) {
7262972d61bSStefano Zampini             for (k=j+1;k<B_N;k++) {
7272972d61bSStefano Zampini               array[k*B_N+j+cum] = array[j*B_N+k+cum];
7282972d61bSStefano Zampini             }
7292972d61bSStefano Zampini           }
7302972d61bSStefano Zampini         }
7319552c7c7SStefano Zampini         if (invert_Stildas) { /* Stildas can be singular */
7322972d61bSStefano Zampini           if (!sub_schurs->is_hermitian) {
73308122e43SStefano Zampini             PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,array2+cum,&B_N,pivots,&B_ierr));
73408122e43SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
73506a4b1faSStefano Zampini             PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,array2+cum,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
73608122e43SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
7372972d61bSStefano Zampini           } else {
7382972d61bSStefano Zampini             PetscInt j,k;
7392972d61bSStefano Zampini 
7402972d61bSStefano Zampini             PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,array2+cum,&B_N,&B_ierr));
7412972d61bSStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
7422972d61bSStefano Zampini             PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,array2+cum,&B_N,&B_ierr));
7432972d61bSStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
7442972d61bSStefano Zampini             for (j=0;j<B_N;j++) {
7452972d61bSStefano Zampini               for (k=j+1;k<B_N;k++) {
7462972d61bSStefano Zampini                 array2[k*B_N+j+cum] = array2[j*B_N+k+cum];
7472972d61bSStefano Zampini               }
7482972d61bSStefano Zampini             }
7492972d61bSStefano Zampini           }
7509552c7c7SStefano Zampini         }
75108122e43SStefano Zampini         ierr = PetscFPTrapPop();CHKERRQ(ierr);
75208122e43SStefano Zampini       }
75308122e43SStefano Zampini       cum += subset_size*subset_size;
75408122e43SStefano Zampini     }
75508122e43SStefano Zampini     ierr = MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_inv_all,&array);CHKERRQ(ierr);
75608122e43SStefano Zampini     ierr = MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_tilda_all,&array2);CHKERRQ(ierr);
75708122e43SStefano Zampini   }
7583202ece2SStefano Zampini 
75906a4b1faSStefano Zampini   ierr = PetscFree2(Bwork,pivots);CHKERRQ(ierr);
760a1337663SStefano Zampini   ierr = MatDestroy(&global_schur_subsets);CHKERRQ(ierr);
761a1337663SStefano Zampini   ierr = MatDestroy(&S_Ej_tilda_all);CHKERRQ(ierr);
76208122e43SStefano Zampini   ierr = MatDestroy(&S_Ej_inv_all);CHKERRQ(ierr);
7633202ece2SStefano Zampini   ierr = MatDestroy(&work_mat);CHKERRQ(ierr);
7645db18549SStefano Zampini   ierr = ISDestroy(&temp_is);CHKERRQ(ierr);
765b1b3d7a2SStefano Zampini   PetscFunctionReturn(0);
766b1b3d7a2SStefano Zampini }
767b1b3d7a2SStefano Zampini 
768b1b3d7a2SStefano Zampini #undef __FUNCT__
769b1b3d7a2SStefano Zampini #define __FUNCT__ "PCBDDCSubSchursInit"
7705db18549SStefano Zampini PetscErrorCode PCBDDCSubSchursInit(PCBDDCSubSchurs sub_schurs, Mat A, Mat S, IS is_I, IS is_B, PCBDDCGraph graph, ISLocalToGlobalMapping BtoNmap, PetscInt seqthreshold)
771b1b3d7a2SStefano Zampini {
7729bb4a8caSStefano Zampini   IS                  *faces,*edges,*all_cc,vertices;
773b1b3d7a2SStefano Zampini   PetscInt            *index_sequential,*index_parallel;
774b1b3d7a2SStefano Zampini   PetscInt            *auxlocal_sequential,*auxlocal_parallel;
775b1b3d7a2SStefano Zampini   PetscInt            *auxglobal_sequential,*auxglobal_parallel;
776b96c3477SStefano Zampini   PetscInt            *auxmapping;
777b1b3d7a2SStefano Zampini   PetscInt            i,max_subset_size;
778b1b3d7a2SStefano Zampini   PetscInt            n_sequential_problems,n_local_sequential_problems,n_parallel_problems,n_local_parallel_problems;
779b1b3d7a2SStefano Zampini   PetscInt            n_faces,n_edges,n_all_cc;
780b1b3d7a2SStefano Zampini   PetscBool           is_sorted;
781b1b3d7a2SStefano Zampini   PetscErrorCode  ierr;
782b1b3d7a2SStefano Zampini 
783b1b3d7a2SStefano Zampini   PetscFunctionBegin;
784b1b3d7a2SStefano Zampini   ierr = ISSorted(is_I,&is_sorted);CHKERRQ(ierr);
785b1b3d7a2SStefano Zampini   if (!is_sorted) {
786b1b3d7a2SStefano Zampini     SETERRQ(PetscObjectComm((PetscObject)is_I),PETSC_ERR_PLIB,"IS for I dofs should be shorted");
787b1b3d7a2SStefano Zampini   }
788b1b3d7a2SStefano Zampini   ierr = ISSorted(is_B,&is_sorted);CHKERRQ(ierr);
789b1b3d7a2SStefano Zampini   if (!is_sorted) {
790b1b3d7a2SStefano Zampini     SETERRQ(PetscObjectComm((PetscObject)is_B),PETSC_ERR_PLIB,"IS for B dofs should be shorted");
791b1b3d7a2SStefano Zampini   }
792b1b3d7a2SStefano Zampini 
793b1b3d7a2SStefano Zampini   /* reset any previous data */
794b1b3d7a2SStefano Zampini   ierr = PCBDDCSubSchursReset(sub_schurs);CHKERRQ(ierr);
795b1b3d7a2SStefano Zampini 
796b1b3d7a2SStefano Zampini   /* get index sets for faces and edges */
7979bb4a8caSStefano Zampini   ierr = PCBDDCGraphGetCandidatesIS(graph,&n_faces,&faces,&n_edges,&edges,&vertices);CHKERRQ(ierr);
798b1b3d7a2SStefano Zampini   n_all_cc = n_faces+n_edges;
79908122e43SStefano Zampini   ierr = PetscBTCreate(n_all_cc,&sub_schurs->is_edge);CHKERRQ(ierr);
80008122e43SStefano Zampini   ierr = PetscBTCreate(n_all_cc,&sub_schurs->computed_Stilda_subs);CHKERRQ(ierr);
801b1b3d7a2SStefano Zampini   ierr = PetscMalloc1(n_all_cc,&all_cc);CHKERRQ(ierr);
802b1b3d7a2SStefano Zampini   for (i=0;i<n_faces;i++) {
803b1b3d7a2SStefano Zampini     all_cc[i] = faces[i];
804b1b3d7a2SStefano Zampini   }
805b1b3d7a2SStefano Zampini   for (i=0;i<n_edges;i++) {
806b1b3d7a2SStefano Zampini     all_cc[n_faces+i] = edges[i];
80708122e43SStefano Zampini     ierr = PetscBTSet(sub_schurs->is_edge,n_faces+i);CHKERRQ(ierr);
808b1b3d7a2SStefano Zampini   }
809b1b3d7a2SStefano Zampini   ierr = PetscFree(faces);CHKERRQ(ierr);
810b1b3d7a2SStefano Zampini   ierr = PetscFree(edges);CHKERRQ(ierr);
811b1b3d7a2SStefano Zampini 
812b1b3d7a2SStefano Zampini   /* map interface's subsets */
813b1b3d7a2SStefano Zampini   max_subset_size = 0;
814b1b3d7a2SStefano Zampini   for (i=0;i<n_all_cc;i++) {
815b1b3d7a2SStefano Zampini     PetscInt subset_size;
816b1b3d7a2SStefano Zampini     ierr = ISGetLocalSize(all_cc[i],&subset_size);CHKERRQ(ierr);
817b1b3d7a2SStefano Zampini     max_subset_size = PetscMax(max_subset_size,subset_size);
818b1b3d7a2SStefano Zampini   }
819b1b3d7a2SStefano Zampini   ierr = PetscMalloc1(max_subset_size,&auxmapping);CHKERRQ(ierr);
820b1b3d7a2SStefano Zampini   ierr = PetscMalloc2(graph->ncc,&auxlocal_sequential,
821b1b3d7a2SStefano Zampini                       graph->ncc,&auxlocal_parallel);CHKERRQ(ierr);
822b1b3d7a2SStefano Zampini   ierr = PetscMalloc2(graph->ncc,&index_sequential,
823b1b3d7a2SStefano Zampini                       graph->ncc,&index_parallel);CHKERRQ(ierr);
824b1b3d7a2SStefano Zampini 
825883469d8SStefano Zampini   /* if threshold is negative uses all sequential problems (possibly using MUMPS) */
826883469d8SStefano Zampini   sub_schurs->use_mumps = PETSC_FALSE;
827883469d8SStefano Zampini   if (seqthreshold < 0) {
828883469d8SStefano Zampini     seqthreshold = max_subset_size;
829883469d8SStefano Zampini #if defined(PETSC_HAVE_MUMPS)
830*4c6709b3SStefano Zampini     sub_schurs->use_mumps = (PetscBool)(!!A);
831883469d8SStefano Zampini #endif
832883469d8SStefano Zampini   }
833883469d8SStefano Zampini 
834b1b3d7a2SStefano Zampini 
835b1b3d7a2SStefano Zampini   /* determine which problem has to be solved in parallel or sequentially */
836b1b3d7a2SStefano Zampini   n_local_sequential_problems = 0;
837b1b3d7a2SStefano Zampini   n_local_parallel_problems = 0;
838b1b3d7a2SStefano Zampini   for (i=0;i<n_all_cc;i++) {
839b1b3d7a2SStefano Zampini     PetscInt       subset_size,j,min_loc = 0;
840b1b3d7a2SStefano Zampini     const PetscInt *idxs;
841b1b3d7a2SStefano Zampini 
842b1b3d7a2SStefano Zampini     ierr = ISGetLocalSize(all_cc[i],&subset_size);CHKERRQ(ierr);
843b1b3d7a2SStefano Zampini     ierr = ISGetIndices(all_cc[i],&idxs);CHKERRQ(ierr);
844b1b3d7a2SStefano Zampini     ierr = ISLocalToGlobalMappingApply(graph->l2gmap,subset_size,idxs,auxmapping);CHKERRQ(ierr);
845b1b3d7a2SStefano Zampini     for (j=1;j<subset_size;j++) {
846b1b3d7a2SStefano Zampini       if (auxmapping[j]<auxmapping[min_loc]) {
847b1b3d7a2SStefano Zampini         min_loc = j;
848b1b3d7a2SStefano Zampini       }
849b1b3d7a2SStefano Zampini     }
850b1b3d7a2SStefano Zampini     if (subset_size > seqthreshold) {
851b1b3d7a2SStefano Zampini       index_parallel[n_local_parallel_problems] = i;
852b1b3d7a2SStefano Zampini       auxlocal_parallel[n_local_parallel_problems] = idxs[min_loc];
853b1b3d7a2SStefano Zampini       n_local_parallel_problems++;
854b1b3d7a2SStefano Zampini     } else {
855b1b3d7a2SStefano Zampini       index_sequential[n_local_sequential_problems] = i;
856b1b3d7a2SStefano Zampini       auxlocal_sequential[n_local_sequential_problems] = idxs[min_loc];
857b1b3d7a2SStefano Zampini       n_local_sequential_problems++;
858b1b3d7a2SStefano Zampini     }
859b1b3d7a2SStefano Zampini     ierr = ISRestoreIndices(all_cc[i],&idxs);CHKERRQ(ierr);
860b1b3d7a2SStefano Zampini   }
861b1b3d7a2SStefano Zampini 
862b1b3d7a2SStefano Zampini   /* Number parallel problems */
863b1b3d7a2SStefano Zampini   auxglobal_parallel = 0;
864b1b3d7a2SStefano Zampini   ierr = PCBDDCSubsetNumbering(PetscObjectComm((PetscObject)graph->l2gmap),graph->l2gmap,n_local_parallel_problems,auxlocal_parallel,PETSC_NULL,&n_parallel_problems,&auxglobal_parallel);CHKERRQ(ierr);
865b1b3d7a2SStefano Zampini 
866b1b3d7a2SStefano Zampini   /* Number sequential problems */
867b1b3d7a2SStefano Zampini   auxglobal_sequential = 0;
868b1b3d7a2SStefano Zampini   ierr = PCBDDCSubsetNumbering(PetscObjectComm((PetscObject)graph->l2gmap),graph->l2gmap,n_local_sequential_problems,auxlocal_sequential,PETSC_NULL,&n_sequential_problems,&auxglobal_sequential);CHKERRQ(ierr);
869b1b3d7a2SStefano Zampini 
870b1b3d7a2SStefano Zampini   /* update info in sub_schurs */
871aa83b6aeSStefano Zampini   if (A) {
8729ab7bb16SStefano Zampini     PetscBool isseqaij;
8739ab7bb16SStefano Zampini 
8749ab7bb16SStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
8759ab7bb16SStefano Zampini     if (isseqaij) {
8761e9c79c2SStefano Zampini       ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
8771e9c79c2SStefano Zampini       sub_schurs->A = A;
8789ab7bb16SStefano Zampini     } else { /* SeqBAIJ matrices does not support symmetry checking, SeqSBAIJ does not support MatPermute */
8799ab7bb16SStefano Zampini       ierr = MatConvert(A,MATSEQAIJ,MAT_INITIAL_MATRIX,&sub_schurs->A);CHKERRQ(ierr);
8809ab7bb16SStefano Zampini     }
8811e9c79c2SStefano Zampini   }
882b1b3d7a2SStefano Zampini   ierr = PetscObjectReference((PetscObject)S);CHKERRQ(ierr);
883b1b3d7a2SStefano Zampini   sub_schurs->S = S;
884b1b3d7a2SStefano Zampini   ierr = PetscObjectReference((PetscObject)is_I);CHKERRQ(ierr);
885b1b3d7a2SStefano Zampini   sub_schurs->is_I = is_I;
886b1b3d7a2SStefano Zampini   ierr = PetscObjectReference((PetscObject)is_B);CHKERRQ(ierr);
887b1b3d7a2SStefano Zampini   sub_schurs->is_B = is_B;
8885db18549SStefano Zampini   ierr = PetscObjectReference((PetscObject)graph->l2gmap);CHKERRQ(ierr);
8895db18549SStefano Zampini   sub_schurs->l2gmap = graph->l2gmap;
8905db18549SStefano Zampini   ierr = PetscObjectReference((PetscObject)BtoNmap);CHKERRQ(ierr);
8915db18549SStefano Zampini   sub_schurs->BtoNmap = BtoNmap;
892b1b3d7a2SStefano Zampini   sub_schurs->n_subs_seq = n_local_sequential_problems;
893b1b3d7a2SStefano Zampini   sub_schurs->n_subs_par = n_local_parallel_problems;
894b1b3d7a2SStefano Zampini   sub_schurs->n_subs_seq_g = n_sequential_problems;
895b1b3d7a2SStefano Zampini   sub_schurs->n_subs_par_g = n_parallel_problems;
896b1b3d7a2SStefano Zampini   sub_schurs->n_subs = sub_schurs->n_subs_seq + sub_schurs->n_subs_par;
897b1b3d7a2SStefano Zampini   sub_schurs->is_subs = all_cc;
8989bb4a8caSStefano Zampini   if (!sub_schurs->use_mumps) { /* for adaptive selection */
899b96c3477SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
900b96c3477SStefano Zampini       ierr = ISSort(sub_schurs->is_subs[i]);CHKERRQ(ierr);
901b96c3477SStefano Zampini     }
9029bb4a8caSStefano Zampini   }
9039bb4a8caSStefano Zampini   sub_schurs->is_Ej_com = vertices;
904b1b3d7a2SStefano Zampini   sub_schurs->index_sequential = index_sequential;
905b1b3d7a2SStefano Zampini   sub_schurs->index_parallel = index_parallel;
906b1b3d7a2SStefano Zampini   sub_schurs->auxglobal_sequential = auxglobal_sequential;
907b1b3d7a2SStefano Zampini   sub_schurs->auxglobal_parallel = auxglobal_parallel;
908b1b3d7a2SStefano Zampini 
909b96c3477SStefano Zampini 
910b96c3477SStefano Zampini   /* allocate space for schur complements */
911b76a6013SStefano Zampini   ierr = PetscCalloc1(sub_schurs->n_subs,&sub_schurs->S_Ej);CHKERRQ(ierr);
912b96c3477SStefano Zampini   sub_schurs->S_Ej_all = NULL;
913b96c3477SStefano Zampini   sub_schurs->sum_S_Ej_all = NULL;
91408122e43SStefano Zampini   sub_schurs->sum_S_Ej_inv_all = NULL;
915b96c3477SStefano Zampini   sub_schurs->sum_S_Ej_tilda_all = NULL;
916b96c3477SStefano Zampini   sub_schurs->is_Ej_all = NULL;
917b96c3477SStefano Zampini 
918b1b3d7a2SStefano Zampini   /* free workspace */
919b1b3d7a2SStefano Zampini   ierr = PetscFree(auxmapping);CHKERRQ(ierr);
920b1b3d7a2SStefano Zampini   ierr = PetscFree2(auxlocal_sequential,auxlocal_parallel);CHKERRQ(ierr);
921b1b3d7a2SStefano Zampini 
922b1b3d7a2SStefano Zampini   PetscFunctionReturn(0);
923b1b3d7a2SStefano Zampini }
924b1b3d7a2SStefano Zampini 
925b1b3d7a2SStefano Zampini #undef __FUNCT__
92634a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCSubSchursCreate"
92734a97f8cSStefano Zampini PetscErrorCode PCBDDCSubSchursCreate(PCBDDCSubSchurs *sub_schurs)
92834a97f8cSStefano Zampini {
92934a97f8cSStefano Zampini   PCBDDCSubSchurs schurs_ctx;
93034a97f8cSStefano Zampini   PetscErrorCode  ierr;
93134a97f8cSStefano Zampini 
93234a97f8cSStefano Zampini   PetscFunctionBegin;
93334a97f8cSStefano Zampini   ierr = PetscNew(&schurs_ctx);CHKERRQ(ierr);
9345ff63025SStefano Zampini   schurs_ctx->n_subs = 0;
93534a97f8cSStefano Zampini   *sub_schurs = schurs_ctx;
93634a97f8cSStefano Zampini   PetscFunctionReturn(0);
93734a97f8cSStefano Zampini }
93834a97f8cSStefano Zampini 
93934a97f8cSStefano Zampini #undef __FUNCT__
94034a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCSubSchursDestroy"
94134a97f8cSStefano Zampini PetscErrorCode PCBDDCSubSchursDestroy(PCBDDCSubSchurs *sub_schurs)
94234a97f8cSStefano Zampini {
94334a97f8cSStefano Zampini   PetscErrorCode ierr;
94434a97f8cSStefano Zampini 
94534a97f8cSStefano Zampini   PetscFunctionBegin;
94634a97f8cSStefano Zampini   ierr = PCBDDCSubSchursReset(*sub_schurs);CHKERRQ(ierr);
94734a97f8cSStefano Zampini   ierr = PetscFree(*sub_schurs);CHKERRQ(ierr);
94834a97f8cSStefano Zampini   PetscFunctionReturn(0);
94934a97f8cSStefano Zampini }
95034a97f8cSStefano Zampini 
95134a97f8cSStefano Zampini #undef __FUNCT__
95234a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCSubSchursReset"
95334a97f8cSStefano Zampini PetscErrorCode PCBDDCSubSchursReset(PCBDDCSubSchurs sub_schurs)
95434a97f8cSStefano Zampini {
95534a97f8cSStefano Zampini   PetscInt       i;
95634a97f8cSStefano Zampini   PetscErrorCode ierr;
95734a97f8cSStefano Zampini 
95834a97f8cSStefano Zampini   PetscFunctionBegin;
9591e9c79c2SStefano Zampini   ierr = MatDestroy(&sub_schurs->A);CHKERRQ(ierr);
960b1b3d7a2SStefano Zampini   ierr = MatDestroy(&sub_schurs->S);CHKERRQ(ierr);
961b1b3d7a2SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_I);CHKERRQ(ierr);
962b1b3d7a2SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_B);CHKERRQ(ierr);
9635db18549SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->l2gmap);CHKERRQ(ierr);
9645db18549SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->BtoNmap);CHKERRQ(ierr);
96541c3ba1bSStefano Zampini   ierr = MatDestroy(&sub_schurs->S_Ej_all);CHKERRQ(ierr);
96641c3ba1bSStefano Zampini   ierr = MatDestroy(&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
96708122e43SStefano Zampini   ierr = MatDestroy(&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
968a1337663SStefano Zampini   ierr = MatDestroy(&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
9695db18549SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_Ej_all);CHKERRQ(ierr);
9709bb4a8caSStefano Zampini   ierr = ISDestroy(&sub_schurs->is_Ej_com);CHKERRQ(ierr);
97108122e43SStefano Zampini   ierr = PetscBTDestroy(&sub_schurs->is_edge);CHKERRQ(ierr);
97208122e43SStefano Zampini   ierr = PetscBTDestroy(&sub_schurs->computed_Stilda_subs);CHKERRQ(ierr);
97334a97f8cSStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
974b1b3d7a2SStefano Zampini     ierr = ISDestroy(&sub_schurs->is_subs[i]);CHKERRQ(ierr);
97534a97f8cSStefano Zampini     ierr = MatDestroy(&sub_schurs->S_Ej[i]);CHKERRQ(ierr);
97634a97f8cSStefano Zampini   }
97768270318SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_I_layer);CHKERRQ(ierr);
9785ff63025SStefano Zampini   if (sub_schurs->n_subs) {
979b1b3d7a2SStefano Zampini     ierr = PetscFree(sub_schurs->is_subs);CHKERRQ(ierr);
980b76a6013SStefano Zampini     ierr = PetscFree(sub_schurs->S_Ej);CHKERRQ(ierr);
9813dc780c3SStefano Zampini   }
982b1b3d7a2SStefano Zampini   ierr = PetscFree2(sub_schurs->index_sequential,sub_schurs->index_parallel);CHKERRQ(ierr);
983b1b3d7a2SStefano Zampini   ierr = PetscFree(sub_schurs->auxglobal_sequential);CHKERRQ(ierr);
984b1b3d7a2SStefano Zampini   ierr = PetscFree(sub_schurs->auxglobal_parallel);CHKERRQ(ierr);
98534a97f8cSStefano Zampini   sub_schurs->n_subs = 0;
98634a97f8cSStefano Zampini   PetscFunctionReturn(0);
98734a97f8cSStefano Zampini }
98834a97f8cSStefano Zampini 
98934a97f8cSStefano Zampini #undef __FUNCT__
99034a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCAdjGetNextLayer_Private"
9912a155e38SStefano Zampini PETSC_STATIC_INLINE PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt* queue_tip,PetscInt n_prev,PetscBT touched,PetscInt* xadj,PetscInt* adjncy,PetscInt* n_added)
99234a97f8cSStefano Zampini {
99334a97f8cSStefano Zampini   PetscInt       i,j,n;
99434a97f8cSStefano Zampini   PetscErrorCode ierr;
99534a97f8cSStefano Zampini 
99634a97f8cSStefano Zampini   PetscFunctionBegin;
99734a97f8cSStefano Zampini   n = 0;
99834a97f8cSStefano Zampini   for (i=-n_prev;i<0;i++) {
99934a97f8cSStefano Zampini     PetscInt start_dof = queue_tip[i];
100034a97f8cSStefano Zampini     for (j=xadj[start_dof];j<xadj[start_dof+1];j++) {
100134a97f8cSStefano Zampini       PetscInt dof = adjncy[j];
100234a97f8cSStefano Zampini       if (!PetscBTLookup(touched,dof)) {
100334a97f8cSStefano Zampini         ierr = PetscBTSet(touched,dof);CHKERRQ(ierr);
100434a97f8cSStefano Zampini         queue_tip[n] = dof;
100534a97f8cSStefano Zampini         n++;
100634a97f8cSStefano Zampini       }
100734a97f8cSStefano Zampini     }
100834a97f8cSStefano Zampini   }
100934a97f8cSStefano Zampini   *n_added = n;
101034a97f8cSStefano Zampini   PetscFunctionReturn(0);
101134a97f8cSStefano Zampini }
1012