xref: /petsc/src/ksp/pc/impls/bddc/bddcschurs.c (revision 5ec10c6ac20a3057719a2d68a0b1ec85fc87085b)
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*);
6*5ec10c6aSStefano Zampini static PetscErrorCode PCBDDCComputeExplicitSchur(Mat,PetscBool,MatReuse,Mat*);
73202ece2SStefano Zampini 
83202ece2SStefano Zampini #undef __FUNCT__
93202ece2SStefano Zampini #define __FUNCT__ "PCBDDCComputeExplicitSchur"
10*5ec10c6aSStefano 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;
173202ece2SStefano Zampini   PetscMPIInt    size;
183202ece2SStefano Zampini   PetscErrorCode ierr;
193202ece2SStefano Zampini 
203202ece2SStefano Zampini   PetscFunctionBegin;
213202ece2SStefano Zampini   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)M),&size);CHKERRQ(ierr);
223202ece2SStefano Zampini   if (size != 1) {
233202ece2SStefano Zampini     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not for parallel matrices");
243202ece2SStefano Zampini   }
253202ece2SStefano Zampini   ierr = MatSchurComplementGetSubMatrices(M, NULL, NULL, &B, &C, &D);CHKERRQ(ierr);
263202ece2SStefano Zampini   ierr = MatSchurComplementGetKSP(M, &ksp);CHKERRQ(ierr);
273202ece2SStefano Zampini   ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr);
283202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) pc, PCLU, &isLU);CHKERRQ(ierr);
293202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) pc, PCILU, &isILU);CHKERRQ(ierr);
303202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) pc, PCCHOLESKY, &isCHOL);CHKERRQ(ierr);
313202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) B, MATSEQDENSE, &Bdense);CHKERRQ(ierr);
323202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) C, MATSEQDENSE, &Cdense);CHKERRQ(ierr);
333202ece2SStefano Zampini   if (!Bdense) {
343202ece2SStefano Zampini     ierr = MatConvert(B, MATSEQDENSE, MAT_INITIAL_MATRIX, &Bd);CHKERRQ(ierr);
353202ece2SStefano Zampini   } else {
363202ece2SStefano Zampini     Bd = B;
373202ece2SStefano Zampini   }
383202ece2SStefano Zampini 
393202ece2SStefano Zampini   if (isLU || isILU || isCHOL) {
403202ece2SStefano Zampini     Mat fact;
413202ece2SStefano Zampini 
423202ece2SStefano Zampini     ierr = KSPSetUp(ksp);CHKERRQ(ierr);
433202ece2SStefano Zampini     ierr = PCFactorGetMatrix(pc, &fact);CHKERRQ(ierr);
443202ece2SStefano Zampini     ierr = MatDuplicate(Bd, MAT_DO_NOT_COPY_VALUES, &AinvBd);CHKERRQ(ierr);
453202ece2SStefano Zampini     ierr = MatMatSolve(fact, Bd, AinvBd);CHKERRQ(ierr);
463202ece2SStefano Zampini   } else {
473202ece2SStefano Zampini     Mat Ainvd;
483202ece2SStefano Zampini 
493202ece2SStefano Zampini     ierr = PCComputeExplicitOperator(pc, &Ainvd);CHKERRQ(ierr);
503202ece2SStefano Zampini     ierr = MatMatMult(Ainvd, Bd, MAT_INITIAL_MATRIX, fill, &AinvBd);CHKERRQ(ierr);
513202ece2SStefano Zampini     ierr = MatDestroy(&Ainvd);CHKERRQ(ierr);
523202ece2SStefano Zampini   }
53*5ec10c6aSStefano Zampini   if (!Bdense & !issym) {
543202ece2SStefano Zampini     ierr = MatDestroy(&Bd);CHKERRQ(ierr);
553202ece2SStefano Zampini   }
56*5ec10c6aSStefano Zampini 
57*5ec10c6aSStefano Zampini   if (!issym) {
583202ece2SStefano Zampini     if (!Cdense) {
593202ece2SStefano Zampini       ierr = MatConvert(C, MATSEQDENSE, MAT_INITIAL_MATRIX, &Cd);CHKERRQ(ierr);
603202ece2SStefano Zampini     } else {
613202ece2SStefano Zampini       Cd = C;
623202ece2SStefano Zampini     }
63*5ec10c6aSStefano Zampini     ierr = MatMatMult(Cd, AinvBd, reuse, fill, S);CHKERRQ(ierr);
643202ece2SStefano Zampini     if (!Cdense) {
653202ece2SStefano Zampini       ierr = MatDestroy(&Cd);CHKERRQ(ierr);
663202ece2SStefano Zampini     }
67*5ec10c6aSStefano Zampini   } else {
68*5ec10c6aSStefano Zampini     ierr = MatTransposeMatMult(Bd, AinvBd, reuse, fill, S);CHKERRQ(ierr);
69*5ec10c6aSStefano Zampini     if (!Bdense) {
70*5ec10c6aSStefano Zampini       ierr = MatDestroy(&Bd);CHKERRQ(ierr);
71*5ec10c6aSStefano Zampini     }
72*5ec10c6aSStefano Zampini   }
73*5ec10c6aSStefano Zampini   ierr = MatDestroy(&AinvBd);CHKERRQ(ierr);
743202ece2SStefano Zampini 
753202ece2SStefano Zampini   if (D) {
763202ece2SStefano Zampini     Mat       Dd;
773202ece2SStefano Zampini     PetscBool Ddense;
783202ece2SStefano Zampini 
793202ece2SStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)D,MATSEQDENSE,&Ddense);CHKERRQ(ierr);
803202ece2SStefano Zampini     if (!Ddense) {
813202ece2SStefano Zampini       ierr = MatConvert(D, MATSEQDENSE, MAT_INITIAL_MATRIX, &Dd);CHKERRQ(ierr);
823202ece2SStefano Zampini     } else {
833202ece2SStefano Zampini       Dd = D;
843202ece2SStefano Zampini     }
853202ece2SStefano Zampini     ierr = MatAYPX(*S,-1.0,Dd,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
863202ece2SStefano Zampini     if (!Ddense) {
873202ece2SStefano Zampini       ierr = MatDestroy(&Dd);CHKERRQ(ierr);
883202ece2SStefano Zampini     }
893202ece2SStefano Zampini   } else {
903202ece2SStefano Zampini     ierr = MatScale(*S,-1.0);CHKERRQ(ierr);
913202ece2SStefano Zampini   }
923202ece2SStefano Zampini   PetscFunctionReturn(0);
933202ece2SStefano Zampini }
9434a97f8cSStefano Zampini 
9534a97f8cSStefano Zampini #undef __FUNCT__
961580ed26SStefano Zampini #define __FUNCT__ "PCBDDCSubSchursSetUp"
979552c7c7SStefano 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)
98b1b3d7a2SStefano Zampini {
99b1b3d7a2SStefano Zampini   Mat                    A_II,A_IB,A_BI,A_BB;
100b1b3d7a2SStefano Zampini   Mat                    AE_II,*AE_IE,*AE_EI,*AE_EE;
101883469d8SStefano Zampini   Mat                    S_all,global_schur_subsets,work_mat;
10208122e43SStefano Zampini   Mat                    S_Ej_tilda_all,S_Ej_inv_all;
1035db18549SStefano Zampini   ISLocalToGlobalMapping l2gmap_subsets;
1045db18549SStefano Zampini   IS                     is_I,*is_subset_B,temp_is;
105883469d8SStefano Zampini   PetscInt               *nnz,*all_local_idx_G,*all_local_idx_B,*all_local_idx_N;
106*5ec10c6aSStefano Zampini   PetscInt               i,subset_size,max_subset_size_seq;
107883469d8SStefano Zampini   PetscInt               extra,local_size,global_size;
10808122e43SStefano Zampini   PetscBLASInt           B_N,B_ierr,B_lwork,*pivots;
10908122e43SStefano Zampini   PetscScalar            *work;
110b1b3d7a2SStefano Zampini   PetscErrorCode         ierr;
111b1b3d7a2SStefano Zampini 
112b1b3d7a2SStefano Zampini   PetscFunctionBegin;
113b1b3d7a2SStefano Zampini   /* get Schur complement matrices */
114883469d8SStefano Zampini   if (!sub_schurs->use_mumps) {
115b7eb3628SStefano Zampini     if (compute_Stilda) {
1163dc780c3SStefano Zampini       SETERRQ(PetscObjectComm((PetscObject)sub_schurs->l2gmap),PETSC_ERR_SUP,"Adaptive selection of constraints requires MUMPS");
117b7eb3628SStefano Zampini     }
118b1b3d7a2SStefano Zampini     ierr = MatSchurComplementGetSubMatrices(sub_schurs->S,&A_II,NULL,&A_IB,&A_BI,&A_BB);CHKERRQ(ierr);
119b1b3d7a2SStefano Zampini     ierr = PetscMalloc4(sub_schurs->n_subs,&is_subset_B,
120b1b3d7a2SStefano Zampini                         sub_schurs->n_subs,&AE_IE,
121b1b3d7a2SStefano Zampini                         sub_schurs->n_subs,&AE_EI,
122b1b3d7a2SStefano Zampini                         sub_schurs->n_subs,&AE_EE);CHKERRQ(ierr);
123a58a30b4SStefano Zampini   } else {
124a58a30b4SStefano Zampini     is_subset_B = NULL;
125a58a30b4SStefano Zampini     AE_IE = NULL;
126a58a30b4SStefano Zampini     AE_EI = NULL;
127a58a30b4SStefano Zampini     AE_EE = NULL;
128b1b3d7a2SStefano Zampini   }
129b1b3d7a2SStefano Zampini 
130b1b3d7a2SStefano Zampini   /* determine interior problems */
131b96c3477SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_I_layer);CHKERRQ(ierr);
1323dc780c3SStefano Zampini   ierr = ISGetLocalSize(sub_schurs->is_I,&i);CHKERRQ(ierr);
1333dc780c3SStefano Zampini   if (nlayers >= 0 && i) { /* Interior problems can be different from the original one */
134b1b3d7a2SStefano Zampini     PetscBT                touched;
135b1b3d7a2SStefano Zampini     const PetscInt*        idx_B;
136b1b3d7a2SStefano Zampini     PetscInt               n_I,n_B,n_local_dofs,n_prev_added,j,layer,*local_numbering;
137b1b3d7a2SStefano Zampini 
1383dc780c3SStefano Zampini     if (xadj == NULL || adjncy == NULL) {
1393dc780c3SStefano Zampini       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cannot request layering without adjacency");
1403dc780c3SStefano Zampini     }
141b1b3d7a2SStefano Zampini     /* get sizes */
142b1b3d7a2SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr);
143b1b3d7a2SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_B,&n_B);CHKERRQ(ierr);
144b1b3d7a2SStefano Zampini 
145b1b3d7a2SStefano Zampini     ierr = PetscMalloc1(n_I+n_B,&local_numbering);CHKERRQ(ierr);
146b1b3d7a2SStefano Zampini     ierr = PetscBTCreate(n_I+n_B,&touched);CHKERRQ(ierr);
147b1b3d7a2SStefano Zampini     ierr = PetscBTMemzero(n_I+n_B,touched);CHKERRQ(ierr);
148b1b3d7a2SStefano Zampini 
149b1b3d7a2SStefano Zampini     /* all boundary dofs must be skipped when adding layers */
150b1b3d7a2SStefano Zampini     ierr = ISGetIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr);
151b1b3d7a2SStefano Zampini     for (j=0;j<n_B;j++) {
152b1b3d7a2SStefano Zampini       ierr = PetscBTSet(touched,idx_B[j]);CHKERRQ(ierr);
153b1b3d7a2SStefano Zampini     }
154b1b3d7a2SStefano Zampini     ierr = PetscMemcpy(local_numbering,idx_B,n_B*sizeof(PetscInt));CHKERRQ(ierr);
155b1b3d7a2SStefano Zampini     ierr = ISRestoreIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr);
156b1b3d7a2SStefano Zampini 
157b1b3d7a2SStefano Zampini     /* add prescribed number of layers of dofs */
158b1b3d7a2SStefano Zampini     n_local_dofs = n_B;
159b1b3d7a2SStefano Zampini     n_prev_added = n_B;
160b1b3d7a2SStefano Zampini     for (layer=0;layer<nlayers;layer++) {
161b1b3d7a2SStefano Zampini       PetscInt n_added;
162b1b3d7a2SStefano Zampini       if (n_local_dofs == n_I+n_B) break;
163b1b3d7a2SStefano Zampini       if (n_local_dofs > n_I+n_B) {
164b1b3d7a2SStefano 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);
165b1b3d7a2SStefano Zampini       }
166b1b3d7a2SStefano Zampini       ierr = PCBDDCAdjGetNextLayer_Private(local_numbering+n_local_dofs,n_prev_added,touched,xadj,adjncy,&n_added);CHKERRQ(ierr);
167b1b3d7a2SStefano Zampini       n_prev_added = n_added;
168b1b3d7a2SStefano Zampini       n_local_dofs += n_added;
169b1b3d7a2SStefano Zampini       if (!n_added) break;
170b1b3d7a2SStefano Zampini     }
171b1b3d7a2SStefano Zampini     ierr = PetscBTDestroy(&touched);CHKERRQ(ierr);
172b1b3d7a2SStefano Zampini 
173883469d8SStefano Zampini     /* IS for I layer dofs in original numbering */
17468270318SStefano 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);
175b1b3d7a2SStefano Zampini     ierr = PetscFree(local_numbering);CHKERRQ(ierr);
17668270318SStefano Zampini     ierr = ISSort(sub_schurs->is_I_layer);CHKERRQ(ierr);
177883469d8SStefano Zampini     /* IS for I layer dofs in I numbering */
178883469d8SStefano Zampini     if (!sub_schurs->use_mumps) {
179b1b3d7a2SStefano Zampini       ISLocalToGlobalMapping ItoNmap;
180b1b3d7a2SStefano Zampini       ierr = ISLocalToGlobalMappingCreateIS(sub_schurs->is_I,&ItoNmap);CHKERRQ(ierr);
18168270318SStefano Zampini       ierr = ISGlobalToLocalMappingApplyIS(ItoNmap,IS_GTOLM_DROP,sub_schurs->is_I_layer,&is_I);CHKERRQ(ierr);
182b1b3d7a2SStefano Zampini       ierr = ISLocalToGlobalMappingDestroy(&ItoNmap);CHKERRQ(ierr);
183b1b3d7a2SStefano Zampini 
184b1b3d7a2SStefano Zampini       /* II block */
185b1b3d7a2SStefano Zampini       ierr = MatGetSubMatrix(A_II,is_I,is_I,MAT_INITIAL_MATRIX,&AE_II);CHKERRQ(ierr);
186b1b3d7a2SStefano Zampini     }
187b1b3d7a2SStefano Zampini   } else {
188b1b3d7a2SStefano Zampini     PetscInt n_I;
189b1b3d7a2SStefano Zampini 
190b1b3d7a2SStefano Zampini     /* IS for I dofs in original numbering */
191b1b3d7a2SStefano Zampini     ierr = PetscObjectReference((PetscObject)sub_schurs->is_I);CHKERRQ(ierr);
19268270318SStefano Zampini     sub_schurs->is_I_layer = sub_schurs->is_I;
193b1b3d7a2SStefano Zampini 
194b1b3d7a2SStefano Zampini     /* IS for I dofs in I numbering (strided 1) */
195883469d8SStefano Zampini     if (!sub_schurs->use_mumps) {
196b1b3d7a2SStefano Zampini       ierr = ISGetSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr);
197b1b3d7a2SStefano Zampini       ierr = ISCreateStride(PetscObjectComm((PetscObject)sub_schurs->is_I),n_I,0,1,&is_I);CHKERRQ(ierr);
198b1b3d7a2SStefano Zampini 
199b1b3d7a2SStefano Zampini       /* II block is the same */
200b1b3d7a2SStefano Zampini       ierr = PetscObjectReference((PetscObject)A_II);CHKERRQ(ierr);
201b1b3d7a2SStefano Zampini       AE_II = A_II;
202b1b3d7a2SStefano Zampini     }
203b1b3d7a2SStefano Zampini   }
204883469d8SStefano Zampini   /* Get info on subset sizes and sum of all subsets sizes */
205*5ec10c6aSStefano Zampini   max_subset_size_seq = 0;
206883469d8SStefano Zampini   local_size = 0;
207883469d8SStefano Zampini   for (i=0;i<sub_schurs->n_subs_seq;i++) {
208883469d8SStefano Zampini     PetscInt j = sub_schurs->index_sequential[i];
209b96c3477SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_subs[j],&subset_size);CHKERRQ(ierr);
210*5ec10c6aSStefano Zampini     max_subset_size_seq = PetscMax(subset_size,max_subset_size_seq);
211883469d8SStefano Zampini     local_size += subset_size;
212883469d8SStefano Zampini   }
213883469d8SStefano Zampini 
214883469d8SStefano Zampini   /* Work arrays for local indices */
215883469d8SStefano Zampini   ierr = PetscMalloc1(local_size,&all_local_idx_B);CHKERRQ(ierr);
216883469d8SStefano Zampini   extra = 0;
217883469d8SStefano Zampini   if (sub_schurs->use_mumps) {
218883469d8SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_I_layer,&extra);CHKERRQ(ierr);
219883469d8SStefano Zampini   }
220883469d8SStefano Zampini   ierr = PetscMalloc1(local_size+extra,&all_local_idx_N);CHKERRQ(ierr);
221883469d8SStefano Zampini   if (extra) {
222883469d8SStefano Zampini     const PetscInt *idxs;
223883469d8SStefano Zampini     ierr = ISGetIndices(sub_schurs->is_I_layer,&idxs);CHKERRQ(ierr);
224883469d8SStefano Zampini     ierr = PetscMemcpy(all_local_idx_N,idxs,extra*sizeof(PetscInt));CHKERRQ(ierr);
225883469d8SStefano Zampini     ierr = ISRestoreIndices(sub_schurs->is_I_layer,&idxs);CHKERRQ(ierr);
226883469d8SStefano Zampini   }
227883469d8SStefano Zampini   ierr = PetscMalloc1(local_size,&nnz);CHKERRQ(ierr);
228883469d8SStefano Zampini 
229883469d8SStefano Zampini   /* Get local indices in local numbering */
230883469d8SStefano Zampini   local_size = 0;
231883469d8SStefano Zampini   for (i=0;i<sub_schurs->n_subs_seq;i++) {
232883469d8SStefano Zampini     PetscInt j;
233883469d8SStefano Zampini     const    PetscInt *idxs;
234883469d8SStefano Zampini 
235883469d8SStefano Zampini     PetscInt local_problem_index = sub_schurs->index_sequential[i];
236b96c3477SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_subs[local_problem_index],&subset_size);CHKERRQ(ierr);
237b96c3477SStefano Zampini     ierr = ISGetIndices(sub_schurs->is_subs[local_problem_index],&idxs);CHKERRQ(ierr);
238883469d8SStefano Zampini     /* subset indices in local numbering */
239883469d8SStefano Zampini     ierr = PetscMemcpy(all_local_idx_N+local_size+extra,idxs,subset_size*sizeof(PetscInt));CHKERRQ(ierr);
240b96c3477SStefano Zampini     ierr = ISRestoreIndices(sub_schurs->is_subs[local_problem_index],&idxs);CHKERRQ(ierr);
241883469d8SStefano Zampini     for (j=0;j<subset_size;j++) nnz[local_size+j] = subset_size;
242883469d8SStefano Zampini     local_size += subset_size;
243883469d8SStefano Zampini   }
244883469d8SStefano Zampini 
245883469d8SStefano Zampini   S_all = NULL;
2469ab7bb16SStefano Zampini   sub_schurs->is_hermitian = PETSC_FALSE;
2479ab7bb16SStefano Zampini   sub_schurs->is_posdef = PETSC_FALSE;
2489552c7c7SStefano Zampini   if (sub_schurs->A) {
2499ab7bb16SStefano Zampini     ierr = MatIsHermitian(sub_schurs->A,0.0,&sub_schurs->is_hermitian);CHKERRQ(ierr);
2509ab7bb16SStefano Zampini     if (sub_schurs->is_hermitian) {
2519ab7bb16SStefano Zampini       PetscScalar val;
2529ab7bb16SStefano Zampini       Vec         vec1,vec2;
2539ab7bb16SStefano Zampini 
2549ab7bb16SStefano Zampini       ierr = MatCreateVecs(sub_schurs->A,&vec1,&vec2);CHKERRQ(ierr);
2559ab7bb16SStefano Zampini       ierr = VecSetRandom(vec1,NULL);
2569ab7bb16SStefano Zampini       ierr = VecCopy(vec1,vec2);CHKERRQ(ierr);
2579ab7bb16SStefano Zampini       ierr = MatMult(sub_schurs->A,vec2,vec1);CHKERRQ(ierr);
2589ab7bb16SStefano Zampini       ierr = VecDot(vec1,vec2,&val);CHKERRQ(ierr);
2599ab7bb16SStefano Zampini       if (PetscRealPart(val) > 0. && PetscImaginaryPart(val) == 0.) sub_schurs->is_posdef = PETSC_TRUE;
2609ab7bb16SStefano Zampini       ierr = VecDestroy(&vec1);CHKERRQ(ierr);
2619ab7bb16SStefano Zampini       ierr = VecDestroy(&vec2);CHKERRQ(ierr);
2629ab7bb16SStefano Zampini     }
2639552c7c7SStefano Zampini   }
2643dc780c3SStefano Zampini   if (sub_schurs->n_subs) {
265883469d8SStefano Zampini     if (!sub_schurs->use_mumps) {
266883469d8SStefano Zampini       /* subsets in original and boundary numbering */
267883469d8SStefano Zampini       for (i=0;i<sub_schurs->n_subs;i++) {
268b96c3477SStefano Zampini         ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,sub_schurs->is_subs[i],&is_subset_B[i]);CHKERRQ(ierr);
269b1b3d7a2SStefano Zampini       }
270b1b3d7a2SStefano Zampini 
271b1b3d7a2SStefano Zampini       /* EE block */
272b1b3d7a2SStefano Zampini       for (i=0;i<sub_schurs->n_subs;i++) {
273b1b3d7a2SStefano Zampini         ierr = MatGetSubMatrix(A_BB,is_subset_B[i],is_subset_B[i],MAT_INITIAL_MATRIX,&AE_EE[i]);CHKERRQ(ierr);
274b1b3d7a2SStefano Zampini       }
275b1b3d7a2SStefano Zampini       /* IE block */
276b1b3d7a2SStefano Zampini       for (i=0;i<sub_schurs->n_subs;i++) {
277b1b3d7a2SStefano Zampini         ierr = MatGetSubMatrix(A_IB,is_I,is_subset_B[i],MAT_INITIAL_MATRIX,&AE_IE[i]);CHKERRQ(ierr);
278b1b3d7a2SStefano Zampini       }
279b1b3d7a2SStefano Zampini       /* EI block */
280b1b3d7a2SStefano Zampini       for (i=0;i<sub_schurs->n_subs;i++) {
281b1b3d7a2SStefano Zampini         ierr = MatGetSubMatrix(A_BI,is_subset_B[i],is_I,MAT_INITIAL_MATRIX,&AE_EI[i]);CHKERRQ(ierr);
282b1b3d7a2SStefano Zampini       }
283b1b3d7a2SStefano Zampini 
284b1b3d7a2SStefano Zampini       /* setup Schur complements on subset */
285b1b3d7a2SStefano Zampini       for (i=0;i<sub_schurs->n_subs;i++) {
286b96c3477SStefano Zampini         ierr = MatDestroy(&sub_schurs->S_Ej[i]);CHKERRQ(ierr);
287b1b3d7a2SStefano Zampini         ierr = MatCreateSchurComplement(AE_II,AE_II,AE_IE[i],AE_EI[i],AE_EE[i],&sub_schurs->S_Ej[i]);CHKERRQ(ierr);
288a58a30b4SStefano Zampini         ierr = MatDestroy(&AE_EE[i]);CHKERRQ(ierr);
289a58a30b4SStefano Zampini         ierr = MatDestroy(&AE_IE[i]);CHKERRQ(ierr);
290a58a30b4SStefano Zampini         ierr = MatDestroy(&AE_EI[i]);CHKERRQ(ierr);
291b1b3d7a2SStefano Zampini         if (AE_II == A_II) { /* we can reuse the same ksp */
292b1b3d7a2SStefano Zampini           KSP ksp;
293b1b3d7a2SStefano Zampini           ierr = MatSchurComplementGetKSP(sub_schurs->S,&ksp);CHKERRQ(ierr);
294b1b3d7a2SStefano Zampini           ierr = MatSchurComplementSetKSP(sub_schurs->S_Ej[i],ksp);CHKERRQ(ierr);
295b1b3d7a2SStefano Zampini         } else { /* build new ksp object which inherits ksp and pc types from the original one */
296b1b3d7a2SStefano Zampini           KSP      origksp,schurksp;
297b1b3d7a2SStefano Zampini           PC       origpc,schurpc;
298b1b3d7a2SStefano Zampini           KSPType  ksp_type;
299b1b3d7a2SStefano Zampini           PCType   pc_type;
300b1b3d7a2SStefano Zampini           PetscInt n_internal;
301b1b3d7a2SStefano Zampini 
302b1b3d7a2SStefano Zampini           ierr = MatSchurComplementGetKSP(sub_schurs->S,&origksp);CHKERRQ(ierr);
303b1b3d7a2SStefano Zampini           ierr = MatSchurComplementGetKSP(sub_schurs->S_Ej[i],&schurksp);CHKERRQ(ierr);
304b1b3d7a2SStefano Zampini           ierr = KSPGetType(origksp,&ksp_type);CHKERRQ(ierr);
305b1b3d7a2SStefano Zampini           ierr = KSPSetType(schurksp,ksp_type);CHKERRQ(ierr);
306b1b3d7a2SStefano Zampini           ierr = KSPGetPC(schurksp,&schurpc);CHKERRQ(ierr);
307b1b3d7a2SStefano Zampini           ierr = KSPGetPC(origksp,&origpc);CHKERRQ(ierr);
308b1b3d7a2SStefano Zampini           ierr = PCGetType(origpc,&pc_type);CHKERRQ(ierr);
309b1b3d7a2SStefano Zampini           ierr = PCSetType(schurpc,pc_type);CHKERRQ(ierr);
310b1b3d7a2SStefano Zampini           ierr = ISGetSize(is_I,&n_internal);CHKERRQ(ierr);
311b1b3d7a2SStefano Zampini           if (n_internal) { /* UMFPACK gives error with 0 sized problems */
312b1b3d7a2SStefano Zampini             MatSolverPackage solver=NULL;
313b1b3d7a2SStefano Zampini             ierr = PCFactorGetMatSolverPackage(origpc,(const MatSolverPackage*)&solver);CHKERRQ(ierr);
314b1b3d7a2SStefano Zampini             if (solver) {
315b1b3d7a2SStefano Zampini               ierr = PCFactorSetMatSolverPackage(schurpc,solver);CHKERRQ(ierr);
316b1b3d7a2SStefano Zampini             }
317b1b3d7a2SStefano Zampini           }
318b1b3d7a2SStefano Zampini           ierr = KSPSetUp(schurksp);CHKERRQ(ierr);
319b1b3d7a2SStefano Zampini         }
320b1b3d7a2SStefano Zampini       }
321b1b3d7a2SStefano Zampini       /* free */
322b1b3d7a2SStefano Zampini       ierr = ISDestroy(&is_I);CHKERRQ(ierr);
323b1b3d7a2SStefano Zampini       ierr = MatDestroy(&AE_II);CHKERRQ(ierr);
324b1b3d7a2SStefano Zampini       for (i=0;i<sub_schurs->n_subs;i++) {
325b1b3d7a2SStefano Zampini         ierr = ISDestroy(&is_subset_B[i]);CHKERRQ(ierr);
326b1b3d7a2SStefano Zampini       }
327b1b3d7a2SStefano Zampini       ierr = PetscFree4(is_subset_B,AE_IE,AE_EI,AE_EE);CHKERRQ(ierr);
328883469d8SStefano Zampini #if defined(PETSC_HAVE_MUMPS)
329883469d8SStefano Zampini     } else {
330883469d8SStefano Zampini       Mat           A,F;
331883469d8SStefano Zampini       IS            is_A_all;
332883469d8SStefano Zampini       PetscInt      *idxs_schur,n_I;
333883469d8SStefano Zampini 
334883469d8SStefano Zampini       /* get working mat */
335883469d8SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_I_layer,&n_I);CHKERRQ(ierr);
336883469d8SStefano Zampini       ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size+n_I,all_local_idx_N,PETSC_COPY_VALUES,&is_A_all);CHKERRQ(ierr);
337b96c3477SStefano Zampini       ierr = MatGetSubMatrixUnsorted(sub_schurs->A,is_A_all,is_A_all,MAT_INITIAL_MATRIX,&A);CHKERRQ(ierr);
338883469d8SStefano Zampini       ierr = ISDestroy(&is_A_all);CHKERRQ(ierr);
339883469d8SStefano Zampini 
34008122e43SStefano Zampini       if (n_I) {
3419ab7bb16SStefano Zampini         if (sub_schurs->is_hermitian && sub_schurs->is_posdef) {
342883469d8SStefano Zampini           ierr = MatGetFactor(A,MATSOLVERMUMPS,MAT_FACTOR_CHOLESKY,&F);CHKERRQ(ierr);
343883469d8SStefano Zampini         } else {
344883469d8SStefano Zampini           ierr = MatGetFactor(A,MATSOLVERMUMPS,MAT_FACTOR_LU,&F);CHKERRQ(ierr);
345883469d8SStefano Zampini         }
346883469d8SStefano Zampini 
347883469d8SStefano Zampini         /* subsets ordered last */
348883469d8SStefano Zampini         ierr = PetscMalloc1(local_size,&idxs_schur);CHKERRQ(ierr);
349883469d8SStefano Zampini         for (i=0;i<local_size;i++) {
350883469d8SStefano Zampini           idxs_schur[i] = n_I+i+1;
351883469d8SStefano Zampini         }
352883469d8SStefano Zampini         ierr = MatMumpsSetSchurIndices(F,local_size,idxs_schur);CHKERRQ(ierr);
353883469d8SStefano Zampini         ierr = PetscFree(idxs_schur);CHKERRQ(ierr);
354883469d8SStefano Zampini 
355883469d8SStefano Zampini         /* factorization step */
3569ab7bb16SStefano Zampini         if (sub_schurs->is_hermitian && sub_schurs->is_posdef) {
357883469d8SStefano Zampini           ierr = MatCholeskyFactorSymbolic(F,A,NULL,NULL);CHKERRQ(ierr);
358883469d8SStefano Zampini           ierr = MatCholeskyFactorNumeric(F,A,NULL);CHKERRQ(ierr);
359883469d8SStefano Zampini         } else {
360883469d8SStefano Zampini           ierr = MatLUFactorSymbolic(F,A,NULL,NULL,NULL);CHKERRQ(ierr);
361883469d8SStefano Zampini           ierr = MatLUFactorNumeric(F,A,NULL);CHKERRQ(ierr);
362883469d8SStefano Zampini         }
363883469d8SStefano Zampini 
364883469d8SStefano Zampini         /* get explicit Schur Complement computed during numeric factorization */
365883469d8SStefano Zampini         ierr = MatMumpsGetSchurComplement(F,&S_all);CHKERRQ(ierr);
366883469d8SStefano Zampini         ierr = MatDestroy(&F);CHKERRQ(ierr);
36708122e43SStefano Zampini       } else {
36808122e43SStefano Zampini         ierr = MatConvert(A,MATSEQDENSE,MAT_INITIAL_MATRIX,&S_all);CHKERRQ(ierr);
36908122e43SStefano Zampini       }
37008122e43SStefano Zampini       ierr = MatDestroy(&A);CHKERRQ(ierr);
371883469d8SStefano Zampini #endif
372b1b3d7a2SStefano Zampini     }
3733dc780c3SStefano Zampini   } else {
3743dc780c3SStefano Zampini     ierr = PetscFree(nnz);CHKERRQ(ierr);
3753dc780c3SStefano Zampini     ierr = PetscFree(all_local_idx_N);CHKERRQ(ierr);
3763dc780c3SStefano Zampini   }
3775db18549SStefano Zampini 
3785db18549SStefano Zampini   if (!sub_schurs->n_subs_seq_g) {
3795db18549SStefano Zampini     PetscFunctionReturn(0);
3805db18549SStefano Zampini   }
3815db18549SStefano Zampini 
3825db18549SStefano Zampini   /* subset indices in local boundary numbering */
383b96c3477SStefano Zampini   if (!sub_schurs->is_Ej_all) {
384883469d8SStefano Zampini     ierr = ISGlobalToLocalMappingApply(sub_schurs->BtoNmap,IS_GTOLM_DROP,local_size,all_local_idx_N+extra,&subset_size,all_local_idx_B);CHKERRQ(ierr);
3855db18549SStefano Zampini     if (subset_size != local_size) {
3865db18549SStefano Zampini       SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in sub_schurs serial (BtoNmap)! %d != %d\n",subset_size,local_size);
3875db18549SStefano Zampini     }
3885db18549SStefano Zampini     ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size,all_local_idx_B,PETSC_OWN_POINTER,&sub_schurs->is_Ej_all);CHKERRQ(ierr);
389b96c3477SStefano Zampini   }
3905db18549SStefano Zampini 
3913202ece2SStefano Zampini   /* Local matrix of all local Schur on subsets transposed */
392b96c3477SStefano Zampini   if (!sub_schurs->S_Ej_all) {
3935db18549SStefano Zampini     ierr = MatCreate(PETSC_COMM_SELF,&sub_schurs->S_Ej_all);CHKERRQ(ierr);
3945db18549SStefano Zampini     ierr = MatSetSizes(sub_schurs->S_Ej_all,PETSC_DECIDE,PETSC_DECIDE,local_size,local_size);CHKERRQ(ierr);
3955db18549SStefano Zampini     ierr = MatSetType(sub_schurs->S_Ej_all,MATAIJ);CHKERRQ(ierr);
3965db18549SStefano Zampini     ierr = MatSeqAIJSetPreallocation(sub_schurs->S_Ej_all,0,nnz);CHKERRQ(ierr);
397b96c3477SStefano Zampini   } else {
398b96c3477SStefano Zampini     ierr = MatZeroEntries(sub_schurs->S_Ej_all);CHKERRQ(ierr);
399b96c3477SStefano Zampini   }
400a1337663SStefano Zampini 
401a1337663SStefano Zampini   S_Ej_tilda_all = 0;
40208122e43SStefano Zampini   S_Ej_inv_all = 0;
403a58a30b4SStefano Zampini   work = NULL;
404a58a30b4SStefano Zampini   pivots = NULL;
40508122e43SStefano Zampini   if (sub_schurs->n_subs && deluxe) { /* workspace needed only for GETRI */
40608122e43SStefano Zampini     PetscScalar lwork;
40708122e43SStefano Zampini 
40808122e43SStefano Zampini     B_lwork = -1;
409*5ec10c6aSStefano Zampini     ierr = PetscBLASIntCast(max_subset_size_seq,&B_N);CHKERRQ(ierr);
41008122e43SStefano Zampini     ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
41108122e43SStefano Zampini     PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,work,&B_N,pivots,&lwork,&B_lwork,&B_ierr));
41208122e43SStefano Zampini     ierr = PetscFPTrapPop();CHKERRQ(ierr);
41308122e43SStefano Zampini     if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to GETRI Lapack routine %d",(int)B_ierr);
41408122e43SStefano Zampini     ierr = PetscBLASIntCast((PetscInt)lwork,&B_lwork);CHKERRQ(ierr);
415*5ec10c6aSStefano Zampini     ierr = PetscMalloc2(B_lwork,&work,max_subset_size_seq,&pivots);CHKERRQ(ierr);
41608122e43SStefano Zampini   }
41708122e43SStefano Zampini 
41808122e43SStefano Zampini   ierr = PetscBTMemzero(sub_schurs->n_subs,sub_schurs->computed_Stilda_subs);CHKERRQ(ierr);
419883469d8SStefano Zampini   if (!sub_schurs->use_mumps) {
420*5ec10c6aSStefano Zampini     if (sub_schurs->n_subs_seq) {
421*5ec10c6aSStefano Zampini       PetscScalar *work;
4225db18549SStefano Zampini       PetscInt    *dummy_idx;
4235db18549SStefano Zampini 
4245db18549SStefano Zampini       /* Work arrays */
425*5ec10c6aSStefano Zampini       ierr = PetscMalloc2(max_subset_size_seq,&dummy_idx,max_subset_size_seq*max_subset_size_seq,&work);CHKERRQ(ierr);
426*5ec10c6aSStefano Zampini 
4273202ece2SStefano Zampini       /* Loop on local problems end compute Schur complements explicitly */
4285db18549SStefano Zampini       local_size = 0;
4295db18549SStefano Zampini       for (i=0;i<sub_schurs->n_subs_seq;i++) {
4303202ece2SStefano Zampini         Mat       S_Ej_expl;
4313202ece2SStefano Zampini         PetscInt  j,lpi = sub_schurs->index_sequential[i];
4323202ece2SStefano Zampini         PetscBool Sdense;
4335db18549SStefano Zampini 
434b96c3477SStefano Zampini         ierr = ISGetLocalSize(sub_schurs->is_subs[lpi],&subset_size);CHKERRQ(ierr);
435*5ec10c6aSStefano Zampini         ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&S_Ej_expl);CHKERRQ(ierr);
436*5ec10c6aSStefano Zampini         /* TODO manage sym flag */
437*5ec10c6aSStefano Zampini         ierr = PCBDDCComputeExplicitSchur(sub_schurs->S_Ej[lpi],PETSC_FALSE,MAT_REUSE_MATRIX,&S_Ej_expl);CHKERRQ(ierr);
4383202ece2SStefano Zampini         ierr = PetscObjectTypeCompare((PetscObject)S_Ej_expl,MATSEQDENSE,&Sdense);CHKERRQ(ierr);
4393202ece2SStefano Zampini         if (Sdense) {
4405db18549SStefano Zampini           for (j=0;j<subset_size;j++) {
4415db18549SStefano Zampini             dummy_idx[j]=local_size+j;
4425db18549SStefano Zampini           }
443*5ec10c6aSStefano Zampini           ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
4443202ece2SStefano Zampini         } else {
4453202ece2SStefano Zampini           SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not yet implemented for sparse matrices");
4465db18549SStefano Zampini         }
4473202ece2SStefano Zampini         local_size += subset_size;
44865d8bf0aSStefano Zampini         ierr = MatDestroy(&sub_schurs->S_Ej[lpi]);CHKERRQ(ierr);
44965d8bf0aSStefano Zampini         sub_schurs->S_Ej[lpi] = S_Ej_expl;
4503202ece2SStefano Zampini       }
451*5ec10c6aSStefano Zampini       ierr = PetscFree2(dummy_idx,work);CHKERRQ(ierr);
452*5ec10c6aSStefano Zampini     }
453883469d8SStefano Zampini   } else {
454883469d8SStefano Zampini     PetscInt    *dummy_idx,n_all;
455*5ec10c6aSStefano Zampini     PetscScalar *work;
456883469d8SStefano Zampini 
457a1337663SStefano Zampini     if (compute_Stilda) {
458a1337663SStefano Zampini       ierr = MatCreate(PETSC_COMM_SELF,&S_Ej_tilda_all);CHKERRQ(ierr);
459a1337663SStefano Zampini       ierr = MatSetSizes(S_Ej_tilda_all,PETSC_DECIDE,PETSC_DECIDE,local_size,local_size);CHKERRQ(ierr);
460a1337663SStefano Zampini       ierr = MatSetType(S_Ej_tilda_all,MATAIJ);CHKERRQ(ierr);
461a1337663SStefano Zampini       ierr = MatSeqAIJSetPreallocation(S_Ej_tilda_all,0,nnz);CHKERRQ(ierr);
46208122e43SStefano Zampini       if (deluxe) {
46308122e43SStefano Zampini         ierr = MatCreate(PETSC_COMM_SELF,&S_Ej_inv_all);CHKERRQ(ierr);
46408122e43SStefano Zampini         ierr = MatSetSizes(S_Ej_inv_all,PETSC_DECIDE,PETSC_DECIDE,local_size,local_size);CHKERRQ(ierr);
46508122e43SStefano Zampini         ierr = MatSetType(S_Ej_inv_all,MATAIJ);CHKERRQ(ierr);
46608122e43SStefano Zampini         ierr = MatSeqAIJSetPreallocation(S_Ej_inv_all,0,nnz);CHKERRQ(ierr);
467a1337663SStefano Zampini       }
46808122e43SStefano Zampini     }
46908122e43SStefano Zampini 
4708275c6a8SStefano Zampini     if (S_all) { /* multilevel */
47108122e43SStefano Zampini       ierr = MatGetSize(S_all,&n_all,NULL);CHKERRQ(ierr);
4728275c6a8SStefano Zampini     } else {
4738275c6a8SStefano Zampini       n_all = 0;
4748275c6a8SStefano Zampini     }
47508122e43SStefano Zampini     local_size = 0;
476883469d8SStefano Zampini     /* Work arrays */
477*5ec10c6aSStefano Zampini     ierr = PetscMalloc2(max_subset_size_seq,&dummy_idx,max_subset_size_seq*max_subset_size_seq,&work);CHKERRQ(ierr);
478883469d8SStefano Zampini 
47965d8bf0aSStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
48065d8bf0aSStefano Zampini       IS  is_E;
48165d8bf0aSStefano Zampini       PetscScalar *vals;
48265d8bf0aSStefano Zampini       PetscInt j;
48365d8bf0aSStefano Zampini 
484b96c3477SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
485a1337663SStefano Zampini       for (j=0;j<subset_size;j++) {
486a1337663SStefano Zampini         dummy_idx[j]=local_size+j;
487a1337663SStefano Zampini       }
48865d8bf0aSStefano Zampini       ierr = ISCreateStride(PetscObjectComm((PetscObject)sub_schurs->is_I),subset_size,local_size,1,&is_E);CHKERRQ(ierr);
489b96c3477SStefano Zampini       if (sub_schurs->S_Ej[i]) {
490b96c3477SStefano Zampini         ierr = MatGetSubMatrix(S_all,is_E,is_E,MAT_REUSE_MATRIX,&sub_schurs->S_Ej[i]);CHKERRQ(ierr);
491b96c3477SStefano Zampini       } else {
492a1337663SStefano Zampini         ierr = MatGetSubMatrix(S_all,is_E,is_E,MAT_INITIAL_MATRIX,&sub_schurs->S_Ej[i]);CHKERRQ(ierr);
493b96c3477SStefano Zampini       }
494a1337663SStefano Zampini       ierr = MatDenseGetArray(sub_schurs->S_Ej[i],&vals);CHKERRQ(ierr);
495a1337663SStefano Zampini       ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,vals,INSERT_VALUES);CHKERRQ(ierr);
496a1337663SStefano Zampini       ierr = MatDenseRestoreArray(sub_schurs->S_Ej[i],&vals);CHKERRQ(ierr);
497a1337663SStefano Zampini 
4986da91a30SStefano Zampini       if (compute_Stilda && ((PetscBTLookup(sub_schurs->is_edge,i) && use_edges) || use_faces)) {
49908122e43SStefano Zampini         ierr = PetscBTSet(sub_schurs->computed_Stilda_subs,i);CHKERRQ(ierr);
50008122e43SStefano Zampini       }
50108122e43SStefano Zampini 
50208122e43SStefano Zampini       if (PetscBTLookup(sub_schurs->computed_Stilda_subs,i)) {
50365d8bf0aSStefano Zampini         Mat Stilda;
50465d8bf0aSStefano Zampini         if (sub_schurs->n_subs == 1) {
50565d8bf0aSStefano Zampini           ierr = PetscObjectReference((PetscObject)sub_schurs->S_Ej[i]);CHKERRQ(ierr);
50665d8bf0aSStefano Zampini           Stilda = sub_schurs->S_Ej[i];
50765d8bf0aSStefano Zampini         } else {
50865d8bf0aSStefano Zampini           KSP ksp;
50965d8bf0aSStefano Zampini           PC  pc;
51065d8bf0aSStefano Zampini           Mat S_EF,S_FE,S_FF,Stilda_impl;
51165d8bf0aSStefano Zampini           IS  is_F;
5122b973097SStefano Zampini           PetscScalar eps=1.e-8;
5132b973097SStefano Zampini           PetscBool chop=PETSC_FALSE;
51465d8bf0aSStefano Zampini 
51565d8bf0aSStefano Zampini           ierr = ISComplement(is_E,0,n_all,&is_F);CHKERRQ(ierr);
51665d8bf0aSStefano Zampini           ierr = MatGetSubMatrix(S_all,is_F,is_F,MAT_INITIAL_MATRIX,&S_FF);CHKERRQ(ierr);
51765d8bf0aSStefano Zampini           ierr = MatGetSubMatrix(S_all,is_F,is_E,MAT_INITIAL_MATRIX,&S_FE);CHKERRQ(ierr);
518*5ec10c6aSStefano Zampini           if (sub_schurs->is_hermitian) { /* just need a placeholder for S_EF */
519*5ec10c6aSStefano Zampini             ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,n_all-subset_size,work,&S_EF);CHKERRQ(ierr);
520*5ec10c6aSStefano Zampini           } else {
521*5ec10c6aSStefano Zampini             ierr = MatGetSubMatrix(S_all,is_E,is_F,MAT_INITIAL_MATRIX,&S_EF);CHKERRQ(ierr);
522*5ec10c6aSStefano Zampini           }
523a1337663SStefano Zampini           ierr = ISDestroy(&is_F);CHKERRQ(ierr);
5242b973097SStefano Zampini           if (chop) {
5252b973097SStefano Zampini             ierr = MatChop(S_FF,eps);CHKERRQ(ierr);
5262b973097SStefano Zampini             ierr = MatConvert(S_FF,MATAIJ,MAT_REUSE_MATRIX,&S_FF);CHKERRQ(ierr);
5272b973097SStefano Zampini           }
528a1337663SStefano Zampini           ierr = MatCreateSchurComplement(S_FF,S_FF,S_FE,S_EF,sub_schurs->S_Ej[i],&Stilda_impl);CHKERRQ(ierr);
529a1337663SStefano Zampini           ierr = MatDestroy(&S_FF);CHKERRQ(ierr);
530a1337663SStefano Zampini           ierr = MatDestroy(&S_FE);CHKERRQ(ierr);
531a1337663SStefano Zampini           ierr = MatDestroy(&S_EF);CHKERRQ(ierr);
53265d8bf0aSStefano Zampini           ierr = MatSchurComplementGetKSP(Stilda_impl,&ksp);CHKERRQ(ierr);
53365d8bf0aSStefano Zampini           ierr = KSPSetType(ksp,KSPPREONLY);CHKERRQ(ierr);
53465d8bf0aSStefano Zampini           ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr);
5359ab7bb16SStefano Zampini           if (sub_schurs->is_hermitian && sub_schurs->is_posdef) {
53665d8bf0aSStefano Zampini             ierr = PCSetType(pc,PCCHOLESKY);CHKERRQ(ierr);
53765d8bf0aSStefano Zampini           } else {
53865d8bf0aSStefano Zampini             ierr = PCSetType(pc,PCLU);CHKERRQ(ierr);
53965d8bf0aSStefano Zampini           }
5402b973097SStefano Zampini           if (!chop) {
54165d8bf0aSStefano Zampini             ierr = PCFactorSetUseInPlace(pc,PETSC_TRUE);CHKERRQ(ierr);
5422b973097SStefano Zampini           } else {
5432b973097SStefano Zampini             ierr = PCFactorSetMatSolverPackage(pc,MATSOLVERMUMPS);CHKERRQ(ierr);
5442b973097SStefano Zampini           }
54565d8bf0aSStefano Zampini           ierr = KSPSetUp(ksp);CHKERRQ(ierr);
546*5ec10c6aSStefano Zampini           ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&Stilda);CHKERRQ(ierr);
547*5ec10c6aSStefano Zampini           ierr = PCBDDCComputeExplicitSchur(Stilda_impl,sub_schurs->is_hermitian,MAT_REUSE_MATRIX,&Stilda);CHKERRQ(ierr);
54865d8bf0aSStefano Zampini           ierr = MatDestroy(&Stilda_impl);CHKERRQ(ierr);
54965d8bf0aSStefano Zampini         }
5502b973097SStefano Zampini /*
5512b973097SStefano Zampini         PetscViewerPushFormat(PETSC_VIEWER_STDOUT_SELF,PETSC_VIEWER_ASCII_MATLAB);
5522b973097SStefano Zampini         ierr = MatView(Stilda,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
5532b973097SStefano Zampini         PetscViewerPopFormat(PETSC_VIEWER_STDOUT_SELF);
5542b973097SStefano Zampini */
555b76a6013SStefano Zampini         ierr = MatDenseGetArray(Stilda,&vals);CHKERRQ(ierr);
55608122e43SStefano Zampini         if (deluxe) { /* when using deluxe scaling, we need (S_1^-1+S_2^-1)^-1 */
55708122e43SStefano Zampini           PetscScalar *vals2;
55808122e43SStefano Zampini 
55908122e43SStefano Zampini           ierr = MatDenseGetArray(sub_schurs->S_Ej[i],&vals2);CHKERRQ(ierr);
56008122e43SStefano Zampini           /* need to be optimized (cholesky) */
56108122e43SStefano Zampini           ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr);
56208122e43SStefano Zampini           ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
5639552c7c7SStefano Zampini           if (invert_Stildas) { /* Stildas can be singular */
56408122e43SStefano Zampini             PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,vals,&B_N,pivots,&B_ierr));
56508122e43SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
56608122e43SStefano Zampini             PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,vals,&B_N,pivots,work,&B_lwork,&B_ierr));
56708122e43SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
5689552c7c7SStefano Zampini           }
56908122e43SStefano Zampini           PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,vals2,&B_N,pivots,&B_ierr));
57008122e43SStefano Zampini           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
57108122e43SStefano Zampini           PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,vals2,&B_N,pivots,work,&B_lwork,&B_ierr));
57208122e43SStefano Zampini           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
57308122e43SStefano Zampini           ierr = PetscFPTrapPop();CHKERRQ(ierr);
57408122e43SStefano Zampini           ierr = MatSetValues(S_Ej_inv_all,subset_size,dummy_idx,subset_size,dummy_idx,vals2,INSERT_VALUES);CHKERRQ(ierr);
57508122e43SStefano Zampini           ierr = MatDenseRestoreArray(sub_schurs->S_Ej[i],&vals2);CHKERRQ(ierr);
57608122e43SStefano Zampini         }
577a1337663SStefano Zampini         ierr = MatSetValues(S_Ej_tilda_all,subset_size,dummy_idx,subset_size,dummy_idx,vals,INSERT_VALUES);CHKERRQ(ierr);
578b76a6013SStefano Zampini         ierr = MatDenseRestoreArray(Stilda,&vals);CHKERRQ(ierr);
579b76a6013SStefano Zampini         ierr = MatDestroy(&Stilda);CHKERRQ(ierr);
58065d8bf0aSStefano Zampini       }
58165d8bf0aSStefano Zampini       ierr = ISDestroy(&is_E);CHKERRQ(ierr);
582883469d8SStefano Zampini       local_size += subset_size;
583883469d8SStefano Zampini     }
584*5ec10c6aSStefano Zampini     ierr = PetscFree2(dummy_idx,work);CHKERRQ(ierr);
5855db18549SStefano Zampini   }
586a1337663SStefano Zampini   ierr = PetscFree(nnz);CHKERRQ(ierr);
587a1337663SStefano Zampini   ierr = MatDestroy(&S_all);CHKERRQ(ierr);
5885db18549SStefano Zampini   ierr = MatAssemblyBegin(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
5895db18549SStefano Zampini   ierr = MatAssemblyEnd(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
5908275c6a8SStefano Zampini   if (S_Ej_tilda_all) {
591a1337663SStefano Zampini     ierr = MatAssemblyBegin(S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
592a1337663SStefano Zampini     ierr = MatAssemblyEnd(S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
5938275c6a8SStefano Zampini   }
5948275c6a8SStefano Zampini   if (S_Ej_inv_all) {
59508122e43SStefano Zampini     ierr = MatAssemblyBegin(S_Ej_inv_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
59608122e43SStefano Zampini     ierr = MatAssemblyEnd(S_Ej_inv_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
59708122e43SStefano Zampini   }
598a1337663SStefano Zampini 
5995db18549SStefano Zampini   /* Global matrix of all assembled Schur on subsets */
600883469d8SStefano 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);
6015db18549SStefano Zampini   ierr = ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)sub_schurs->l2gmap),1,local_size,all_local_idx_G,PETSC_COPY_VALUES,&l2gmap_subsets);CHKERRQ(ierr);
6025db18549SStefano Zampini   ierr = MatCreateIS(PetscObjectComm((PetscObject)sub_schurs->l2gmap),1,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size,l2gmap_subsets,&work_mat);CHKERRQ(ierr);
6035db18549SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&l2gmap_subsets);CHKERRQ(ierr);
6045db18549SStefano Zampini   ierr = MatISSetLocalMat(work_mat,sub_schurs->S_Ej_all);CHKERRQ(ierr);
6053927de2eSStefano Zampini   ierr = MatCreate(PetscObjectComm((PetscObject)work_mat),&global_schur_subsets);CHKERRQ(ierr);
6063927de2eSStefano Zampini   ierr = MatSetSizes(global_schur_subsets,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size);CHKERRQ(ierr);
6073927de2eSStefano Zampini   ierr = MatSetType(global_schur_subsets,MATMPIAIJ);CHKERRQ(ierr);
6083927de2eSStefano Zampini   ierr = MatISSetMPIXAIJPreallocation_Private(work_mat,global_schur_subsets,PETSC_TRUE);CHKERRQ(ierr);
6093927de2eSStefano Zampini   ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
6105db18549SStefano Zampini   /* Get local part of (\sum_j S_Ej) */
611883469d8SStefano Zampini   ierr = PetscFree(all_local_idx_N);CHKERRQ(ierr);
612883469d8SStefano Zampini   ierr = ISCreateGeneral(PetscObjectComm((PetscObject)sub_schurs->l2gmap),local_size,all_local_idx_G,PETSC_OWN_POINTER,&temp_is);CHKERRQ(ierr);
613b96c3477SStefano Zampini   if (sub_schurs->sum_S_Ej_all) {
614b96c3477SStefano Zampini     ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,MAT_REUSE_MATRIX,&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
615b96c3477SStefano Zampini   } else {
616b96c3477SStefano Zampini     ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,MAT_INITIAL_MATRIX,&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
617b96c3477SStefano Zampini   }
61808122e43SStefano Zampini 
6198275c6a8SStefano Zampini   if (S_Ej_tilda_all) {
620a1337663SStefano Zampini     ierr = MatISSetLocalMat(work_mat,S_Ej_tilda_all);CHKERRQ(ierr);
621a1337663SStefano Zampini     ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
622b96c3477SStefano Zampini     if (sub_schurs->sum_S_Ej_tilda_all) {
623b96c3477SStefano Zampini       ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,MAT_REUSE_MATRIX,&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
624b96c3477SStefano Zampini     } else {
625b96c3477SStefano Zampini       ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,MAT_INITIAL_MATRIX,&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
626b96c3477SStefano Zampini     }
6278275c6a8SStefano Zampini   }
6288275c6a8SStefano Zampini   if (S_Ej_inv_all) {
62908122e43SStefano Zampini     PetscInt    cum;
63008122e43SStefano Zampini     PetscScalar *array,*array2;
63108122e43SStefano Zampini     ierr = MatISSetLocalMat(work_mat,S_Ej_inv_all);CHKERRQ(ierr);
63208122e43SStefano Zampini     ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
63308122e43SStefano Zampini     if (sub_schurs->sum_S_Ej_inv_all) {
63408122e43SStefano Zampini       ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,MAT_REUSE_MATRIX,&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
63508122e43SStefano Zampini     } else {
63608122e43SStefano Zampini       ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,MAT_INITIAL_MATRIX,&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
63708122e43SStefano Zampini     }
63808122e43SStefano Zampini     /* invert blocks of sum_S_Ej_inv_all */
63908122e43SStefano Zampini     ierr = MatSeqAIJGetArray(sub_schurs->sum_S_Ej_inv_all,&array);CHKERRQ(ierr);
64008122e43SStefano Zampini     ierr = MatSeqAIJGetArray(sub_schurs->sum_S_Ej_tilda_all,&array2);CHKERRQ(ierr);
64108122e43SStefano Zampini     cum = 0;
64208122e43SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
64308122e43SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
64408122e43SStefano Zampini       if (PetscBTLookup(sub_schurs->computed_Stilda_subs,i)) {
64508122e43SStefano Zampini         /* need to be optimized (cholesky) */
64608122e43SStefano Zampini         ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr);
64708122e43SStefano Zampini         ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
64808122e43SStefano Zampini         PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,array+cum,&B_N,pivots,&B_ierr));
64908122e43SStefano Zampini         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
65008122e43SStefano Zampini         PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,array+cum,&B_N,pivots,work,&B_lwork,&B_ierr));
65108122e43SStefano Zampini         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
6529552c7c7SStefano Zampini         if (invert_Stildas) { /* Stildas can be singular */
65308122e43SStefano Zampini           PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,array2+cum,&B_N,pivots,&B_ierr));
65408122e43SStefano Zampini           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
65508122e43SStefano Zampini           PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,array2+cum,&B_N,pivots,work,&B_lwork,&B_ierr));
65608122e43SStefano Zampini           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
6579552c7c7SStefano Zampini         }
65808122e43SStefano Zampini         ierr = PetscFPTrapPop();CHKERRQ(ierr);
65908122e43SStefano Zampini       }
66008122e43SStefano Zampini       cum += subset_size*subset_size;
66108122e43SStefano Zampini     }
66208122e43SStefano Zampini     ierr = MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_inv_all,&array);CHKERRQ(ierr);
66308122e43SStefano Zampini     ierr = MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_tilda_all,&array2);CHKERRQ(ierr);
66408122e43SStefano Zampini   }
6653202ece2SStefano Zampini 
66608122e43SStefano Zampini   ierr = PetscFree2(work,pivots);CHKERRQ(ierr);
667a1337663SStefano Zampini   ierr = MatDestroy(&global_schur_subsets);CHKERRQ(ierr);
668a1337663SStefano Zampini   ierr = MatDestroy(&S_Ej_tilda_all);CHKERRQ(ierr);
66908122e43SStefano Zampini   ierr = MatDestroy(&S_Ej_inv_all);CHKERRQ(ierr);
6703202ece2SStefano Zampini   ierr = MatDestroy(&work_mat);CHKERRQ(ierr);
6715db18549SStefano Zampini   ierr = ISDestroy(&temp_is);CHKERRQ(ierr);
672b1b3d7a2SStefano Zampini   PetscFunctionReturn(0);
673b1b3d7a2SStefano Zampini }
674b1b3d7a2SStefano Zampini 
675b1b3d7a2SStefano Zampini #undef __FUNCT__
676b1b3d7a2SStefano Zampini #define __FUNCT__ "PCBDDCSubSchursInit"
6775db18549SStefano Zampini PetscErrorCode PCBDDCSubSchursInit(PCBDDCSubSchurs sub_schurs, Mat A, Mat S, IS is_I, IS is_B, PCBDDCGraph graph, ISLocalToGlobalMapping BtoNmap, PetscInt seqthreshold)
678b1b3d7a2SStefano Zampini {
6799bb4a8caSStefano Zampini   IS                  *faces,*edges,*all_cc,vertices;
680b1b3d7a2SStefano Zampini   PetscInt            *index_sequential,*index_parallel;
681b1b3d7a2SStefano Zampini   PetscInt            *auxlocal_sequential,*auxlocal_parallel;
682b1b3d7a2SStefano Zampini   PetscInt            *auxglobal_sequential,*auxglobal_parallel;
683b96c3477SStefano Zampini   PetscInt            *auxmapping;
684b1b3d7a2SStefano Zampini   PetscInt            i,max_subset_size;
685b1b3d7a2SStefano Zampini   PetscInt            n_sequential_problems,n_local_sequential_problems,n_parallel_problems,n_local_parallel_problems;
686b1b3d7a2SStefano Zampini   PetscInt            n_faces,n_edges,n_all_cc;
687b1b3d7a2SStefano Zampini   PetscBool           is_sorted;
688b1b3d7a2SStefano Zampini   PetscErrorCode  ierr;
689b1b3d7a2SStefano Zampini 
690b1b3d7a2SStefano Zampini   PetscFunctionBegin;
691b1b3d7a2SStefano Zampini   ierr = ISSorted(is_I,&is_sorted);CHKERRQ(ierr);
692b1b3d7a2SStefano Zampini   if (!is_sorted) {
693b1b3d7a2SStefano Zampini     SETERRQ(PetscObjectComm((PetscObject)is_I),PETSC_ERR_PLIB,"IS for I dofs should be shorted");
694b1b3d7a2SStefano Zampini   }
695b1b3d7a2SStefano Zampini   ierr = ISSorted(is_B,&is_sorted);CHKERRQ(ierr);
696b1b3d7a2SStefano Zampini   if (!is_sorted) {
697b1b3d7a2SStefano Zampini     SETERRQ(PetscObjectComm((PetscObject)is_B),PETSC_ERR_PLIB,"IS for B dofs should be shorted");
698b1b3d7a2SStefano Zampini   }
699b1b3d7a2SStefano Zampini 
700b1b3d7a2SStefano Zampini   /* reset any previous data */
701b1b3d7a2SStefano Zampini   ierr = PCBDDCSubSchursReset(sub_schurs);CHKERRQ(ierr);
702b1b3d7a2SStefano Zampini 
703b1b3d7a2SStefano Zampini   /* get index sets for faces and edges */
7049bb4a8caSStefano Zampini   ierr = PCBDDCGraphGetCandidatesIS(graph,&n_faces,&faces,&n_edges,&edges,&vertices);CHKERRQ(ierr);
705b1b3d7a2SStefano Zampini   n_all_cc = n_faces+n_edges;
70608122e43SStefano Zampini   ierr = PetscBTCreate(n_all_cc,&sub_schurs->is_edge);CHKERRQ(ierr);
70708122e43SStefano Zampini   ierr = PetscBTCreate(n_all_cc,&sub_schurs->computed_Stilda_subs);CHKERRQ(ierr);
708b1b3d7a2SStefano Zampini   ierr = PetscMalloc1(n_all_cc,&all_cc);CHKERRQ(ierr);
709b1b3d7a2SStefano Zampini   for (i=0;i<n_faces;i++) {
710b1b3d7a2SStefano Zampini     all_cc[i] = faces[i];
711b1b3d7a2SStefano Zampini   }
712b1b3d7a2SStefano Zampini   for (i=0;i<n_edges;i++) {
713b1b3d7a2SStefano Zampini     all_cc[n_faces+i] = edges[i];
71408122e43SStefano Zampini     ierr = PetscBTSet(sub_schurs->is_edge,n_faces+i);CHKERRQ(ierr);
715b1b3d7a2SStefano Zampini   }
716b1b3d7a2SStefano Zampini   ierr = PetscFree(faces);CHKERRQ(ierr);
717b1b3d7a2SStefano Zampini   ierr = PetscFree(edges);CHKERRQ(ierr);
718b1b3d7a2SStefano Zampini 
719b1b3d7a2SStefano Zampini   /* map interface's subsets */
720b1b3d7a2SStefano Zampini   max_subset_size = 0;
721b1b3d7a2SStefano Zampini   for (i=0;i<n_all_cc;i++) {
722b1b3d7a2SStefano Zampini     PetscInt subset_size;
723b1b3d7a2SStefano Zampini     ierr = ISGetLocalSize(all_cc[i],&subset_size);CHKERRQ(ierr);
724b1b3d7a2SStefano Zampini     max_subset_size = PetscMax(max_subset_size,subset_size);
725b1b3d7a2SStefano Zampini   }
726b1b3d7a2SStefano Zampini   ierr = PetscMalloc1(max_subset_size,&auxmapping);CHKERRQ(ierr);
727b1b3d7a2SStefano Zampini   ierr = PetscMalloc2(graph->ncc,&auxlocal_sequential,
728b1b3d7a2SStefano Zampini                       graph->ncc,&auxlocal_parallel);CHKERRQ(ierr);
729b1b3d7a2SStefano Zampini   ierr = PetscMalloc2(graph->ncc,&index_sequential,
730b1b3d7a2SStefano Zampini                       graph->ncc,&index_parallel);CHKERRQ(ierr);
731b1b3d7a2SStefano Zampini 
732883469d8SStefano Zampini   /* if threshold is negative uses all sequential problems (possibly using MUMPS) */
733883469d8SStefano Zampini   sub_schurs->use_mumps = PETSC_FALSE;
734883469d8SStefano Zampini   if (seqthreshold < 0) {
735883469d8SStefano Zampini     seqthreshold = max_subset_size;
736883469d8SStefano Zampini #if defined(PETSC_HAVE_MUMPS)
737883469d8SStefano Zampini     sub_schurs->use_mumps = !!A;
738883469d8SStefano Zampini #endif
739883469d8SStefano Zampini   }
740883469d8SStefano Zampini 
741b1b3d7a2SStefano Zampini 
742b1b3d7a2SStefano Zampini   /* determine which problem has to be solved in parallel or sequentially */
743b1b3d7a2SStefano Zampini   n_local_sequential_problems = 0;
744b1b3d7a2SStefano Zampini   n_local_parallel_problems = 0;
745b1b3d7a2SStefano Zampini   for (i=0;i<n_all_cc;i++) {
746b1b3d7a2SStefano Zampini     PetscInt       subset_size,j,min_loc = 0;
747b1b3d7a2SStefano Zampini     const PetscInt *idxs;
748b1b3d7a2SStefano Zampini 
749b1b3d7a2SStefano Zampini     ierr = ISGetLocalSize(all_cc[i],&subset_size);CHKERRQ(ierr);
750b1b3d7a2SStefano Zampini     ierr = ISGetIndices(all_cc[i],&idxs);CHKERRQ(ierr);
751b1b3d7a2SStefano Zampini     ierr = ISLocalToGlobalMappingApply(graph->l2gmap,subset_size,idxs,auxmapping);CHKERRQ(ierr);
752b1b3d7a2SStefano Zampini     for (j=1;j<subset_size;j++) {
753b1b3d7a2SStefano Zampini       if (auxmapping[j]<auxmapping[min_loc]) {
754b1b3d7a2SStefano Zampini         min_loc = j;
755b1b3d7a2SStefano Zampini       }
756b1b3d7a2SStefano Zampini     }
757b1b3d7a2SStefano Zampini     if (subset_size > seqthreshold) {
758b1b3d7a2SStefano Zampini       index_parallel[n_local_parallel_problems] = i;
759b1b3d7a2SStefano Zampini       auxlocal_parallel[n_local_parallel_problems] = idxs[min_loc];
760b1b3d7a2SStefano Zampini       n_local_parallel_problems++;
761b1b3d7a2SStefano Zampini     } else {
762b1b3d7a2SStefano Zampini       index_sequential[n_local_sequential_problems] = i;
763b1b3d7a2SStefano Zampini       auxlocal_sequential[n_local_sequential_problems] = idxs[min_loc];
764b1b3d7a2SStefano Zampini       n_local_sequential_problems++;
765b1b3d7a2SStefano Zampini     }
766b1b3d7a2SStefano Zampini     ierr = ISRestoreIndices(all_cc[i],&idxs);CHKERRQ(ierr);
767b1b3d7a2SStefano Zampini   }
768b1b3d7a2SStefano Zampini 
769b1b3d7a2SStefano Zampini   /* Number parallel problems */
770b1b3d7a2SStefano Zampini   auxglobal_parallel = 0;
771b1b3d7a2SStefano Zampini   ierr = PCBDDCSubsetNumbering(PetscObjectComm((PetscObject)graph->l2gmap),graph->l2gmap,n_local_parallel_problems,auxlocal_parallel,PETSC_NULL,&n_parallel_problems,&auxglobal_parallel);CHKERRQ(ierr);
772b1b3d7a2SStefano Zampini 
773b1b3d7a2SStefano Zampini   /* Number sequential problems */
774b1b3d7a2SStefano Zampini   auxglobal_sequential = 0;
775b1b3d7a2SStefano Zampini   ierr = PCBDDCSubsetNumbering(PetscObjectComm((PetscObject)graph->l2gmap),graph->l2gmap,n_local_sequential_problems,auxlocal_sequential,PETSC_NULL,&n_sequential_problems,&auxglobal_sequential);CHKERRQ(ierr);
776b1b3d7a2SStefano Zampini 
777b1b3d7a2SStefano Zampini   /* update info in sub_schurs */
778883469d8SStefano Zampini   if (sub_schurs->use_mumps && A) {
7799ab7bb16SStefano Zampini     PetscBool isseqaij;
7809ab7bb16SStefano Zampini 
7819ab7bb16SStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
7829ab7bb16SStefano Zampini     if (isseqaij) {
7831e9c79c2SStefano Zampini       ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
7841e9c79c2SStefano Zampini       sub_schurs->A = A;
7859ab7bb16SStefano Zampini     } else { /* SeqBAIJ matrices does not support symmetry checking, SeqSBAIJ does not support MatPermute */
7869ab7bb16SStefano Zampini       ierr = MatConvert(A,MATSEQAIJ,MAT_INITIAL_MATRIX,&sub_schurs->A);CHKERRQ(ierr);
7879ab7bb16SStefano Zampini     }
7881e9c79c2SStefano Zampini   }
789b1b3d7a2SStefano Zampini   ierr = PetscObjectReference((PetscObject)S);CHKERRQ(ierr);
790b1b3d7a2SStefano Zampini   sub_schurs->S = S;
791b1b3d7a2SStefano Zampini   ierr = PetscObjectReference((PetscObject)is_I);CHKERRQ(ierr);
792b1b3d7a2SStefano Zampini   sub_schurs->is_I = is_I;
793b1b3d7a2SStefano Zampini   ierr = PetscObjectReference((PetscObject)is_B);CHKERRQ(ierr);
794b1b3d7a2SStefano Zampini   sub_schurs->is_B = is_B;
7955db18549SStefano Zampini   ierr = PetscObjectReference((PetscObject)graph->l2gmap);CHKERRQ(ierr);
7965db18549SStefano Zampini   sub_schurs->l2gmap = graph->l2gmap;
7975db18549SStefano Zampini   ierr = PetscObjectReference((PetscObject)BtoNmap);CHKERRQ(ierr);
7985db18549SStefano Zampini   sub_schurs->BtoNmap = BtoNmap;
799b1b3d7a2SStefano Zampini   sub_schurs->n_subs_seq = n_local_sequential_problems;
800b1b3d7a2SStefano Zampini   sub_schurs->n_subs_par = n_local_parallel_problems;
801b1b3d7a2SStefano Zampini   sub_schurs->n_subs_seq_g = n_sequential_problems;
802b1b3d7a2SStefano Zampini   sub_schurs->n_subs_par_g = n_parallel_problems;
803b1b3d7a2SStefano Zampini   sub_schurs->n_subs = sub_schurs->n_subs_seq + sub_schurs->n_subs_par;
804b1b3d7a2SStefano Zampini   sub_schurs->is_subs = all_cc;
8059bb4a8caSStefano Zampini   if (!sub_schurs->use_mumps) { /* for adaptive selection */
806b96c3477SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
807b96c3477SStefano Zampini       ierr = ISSort(sub_schurs->is_subs[i]);CHKERRQ(ierr);
808b96c3477SStefano Zampini     }
8099bb4a8caSStefano Zampini   }
8109bb4a8caSStefano Zampini   sub_schurs->is_Ej_com = vertices;
811b1b3d7a2SStefano Zampini   sub_schurs->index_sequential = index_sequential;
812b1b3d7a2SStefano Zampini   sub_schurs->index_parallel = index_parallel;
813b1b3d7a2SStefano Zampini   sub_schurs->auxglobal_sequential = auxglobal_sequential;
814b1b3d7a2SStefano Zampini   sub_schurs->auxglobal_parallel = auxglobal_parallel;
815b1b3d7a2SStefano Zampini 
816b96c3477SStefano Zampini 
817b96c3477SStefano Zampini   /* allocate space for schur complements */
818b76a6013SStefano Zampini   ierr = PetscCalloc1(sub_schurs->n_subs,&sub_schurs->S_Ej);CHKERRQ(ierr);
819b96c3477SStefano Zampini   sub_schurs->S_Ej_all = NULL;
820b96c3477SStefano Zampini   sub_schurs->sum_S_Ej_all = NULL;
82108122e43SStefano Zampini   sub_schurs->sum_S_Ej_inv_all = NULL;
822b96c3477SStefano Zampini   sub_schurs->sum_S_Ej_tilda_all = NULL;
823b96c3477SStefano Zampini   sub_schurs->is_Ej_all = NULL;
824b96c3477SStefano Zampini 
825b1b3d7a2SStefano Zampini   /* free workspace */
826b1b3d7a2SStefano Zampini   ierr = PetscFree(auxmapping);CHKERRQ(ierr);
827b1b3d7a2SStefano Zampini   ierr = PetscFree2(auxlocal_sequential,auxlocal_parallel);CHKERRQ(ierr);
828b1b3d7a2SStefano Zampini 
829b1b3d7a2SStefano Zampini   PetscFunctionReturn(0);
830b1b3d7a2SStefano Zampini }
831b1b3d7a2SStefano Zampini 
832b1b3d7a2SStefano Zampini #undef __FUNCT__
83334a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCSubSchursCreate"
83434a97f8cSStefano Zampini PetscErrorCode PCBDDCSubSchursCreate(PCBDDCSubSchurs *sub_schurs)
83534a97f8cSStefano Zampini {
83634a97f8cSStefano Zampini   PCBDDCSubSchurs schurs_ctx;
83734a97f8cSStefano Zampini   PetscErrorCode  ierr;
83834a97f8cSStefano Zampini 
83934a97f8cSStefano Zampini   PetscFunctionBegin;
84034a97f8cSStefano Zampini   ierr = PetscNew(&schurs_ctx);CHKERRQ(ierr);
8415ff63025SStefano Zampini   schurs_ctx->n_subs = 0;
84234a97f8cSStefano Zampini   *sub_schurs = schurs_ctx;
84334a97f8cSStefano Zampini   PetscFunctionReturn(0);
84434a97f8cSStefano Zampini }
84534a97f8cSStefano Zampini 
84634a97f8cSStefano Zampini #undef __FUNCT__
84734a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCSubSchursDestroy"
84834a97f8cSStefano Zampini PetscErrorCode PCBDDCSubSchursDestroy(PCBDDCSubSchurs *sub_schurs)
84934a97f8cSStefano Zampini {
85034a97f8cSStefano Zampini   PetscErrorCode ierr;
85134a97f8cSStefano Zampini 
85234a97f8cSStefano Zampini   PetscFunctionBegin;
85334a97f8cSStefano Zampini   ierr = PCBDDCSubSchursReset(*sub_schurs);CHKERRQ(ierr);
85434a97f8cSStefano Zampini   ierr = PetscFree(*sub_schurs);CHKERRQ(ierr);
85534a97f8cSStefano Zampini   PetscFunctionReturn(0);
85634a97f8cSStefano Zampini }
85734a97f8cSStefano Zampini 
85834a97f8cSStefano Zampini #undef __FUNCT__
85934a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCSubSchursReset"
86034a97f8cSStefano Zampini PetscErrorCode PCBDDCSubSchursReset(PCBDDCSubSchurs sub_schurs)
86134a97f8cSStefano Zampini {
86234a97f8cSStefano Zampini   PetscInt       i;
86334a97f8cSStefano Zampini   PetscErrorCode ierr;
86434a97f8cSStefano Zampini 
86534a97f8cSStefano Zampini   PetscFunctionBegin;
8661e9c79c2SStefano Zampini   ierr = MatDestroy(&sub_schurs->A);CHKERRQ(ierr);
867b1b3d7a2SStefano Zampini   ierr = MatDestroy(&sub_schurs->S);CHKERRQ(ierr);
868b1b3d7a2SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_I);CHKERRQ(ierr);
869b1b3d7a2SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_B);CHKERRQ(ierr);
8705db18549SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->l2gmap);CHKERRQ(ierr);
8715db18549SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->BtoNmap);CHKERRQ(ierr);
87241c3ba1bSStefano Zampini   ierr = MatDestroy(&sub_schurs->S_Ej_all);CHKERRQ(ierr);
87341c3ba1bSStefano Zampini   ierr = MatDestroy(&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
87408122e43SStefano Zampini   ierr = MatDestroy(&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
875a1337663SStefano Zampini   ierr = MatDestroy(&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
8765db18549SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_Ej_all);CHKERRQ(ierr);
8779bb4a8caSStefano Zampini   ierr = ISDestroy(&sub_schurs->is_Ej_com);CHKERRQ(ierr);
87808122e43SStefano Zampini   ierr = PetscBTDestroy(&sub_schurs->is_edge);CHKERRQ(ierr);
87908122e43SStefano Zampini   ierr = PetscBTDestroy(&sub_schurs->computed_Stilda_subs);CHKERRQ(ierr);
88034a97f8cSStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
881b1b3d7a2SStefano Zampini     ierr = ISDestroy(&sub_schurs->is_subs[i]);CHKERRQ(ierr);
88234a97f8cSStefano Zampini     ierr = MatDestroy(&sub_schurs->S_Ej[i]);CHKERRQ(ierr);
88334a97f8cSStefano Zampini   }
88468270318SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_I_layer);CHKERRQ(ierr);
8855ff63025SStefano Zampini   if (sub_schurs->n_subs) {
886b1b3d7a2SStefano Zampini     ierr = PetscFree(sub_schurs->is_subs);CHKERRQ(ierr);
887b76a6013SStefano Zampini     ierr = PetscFree(sub_schurs->S_Ej);CHKERRQ(ierr);
8883dc780c3SStefano Zampini   }
889b1b3d7a2SStefano Zampini   ierr = PetscFree2(sub_schurs->index_sequential,sub_schurs->index_parallel);CHKERRQ(ierr);
890b1b3d7a2SStefano Zampini   ierr = PetscFree(sub_schurs->auxglobal_sequential);CHKERRQ(ierr);
891b1b3d7a2SStefano Zampini   ierr = PetscFree(sub_schurs->auxglobal_parallel);CHKERRQ(ierr);
89234a97f8cSStefano Zampini   sub_schurs->n_subs = 0;
89334a97f8cSStefano Zampini   PetscFunctionReturn(0);
89434a97f8cSStefano Zampini }
89534a97f8cSStefano Zampini 
89634a97f8cSStefano Zampini #undef __FUNCT__
89734a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCAdjGetNextLayer_Private"
8982a155e38SStefano Zampini PETSC_STATIC_INLINE PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt* queue_tip,PetscInt n_prev,PetscBT touched,PetscInt* xadj,PetscInt* adjncy,PetscInt* n_added)
89934a97f8cSStefano Zampini {
90034a97f8cSStefano Zampini   PetscInt       i,j,n;
90134a97f8cSStefano Zampini   PetscErrorCode ierr;
90234a97f8cSStefano Zampini 
90334a97f8cSStefano Zampini   PetscFunctionBegin;
90434a97f8cSStefano Zampini   n = 0;
90534a97f8cSStefano Zampini   for (i=-n_prev;i<0;i++) {
90634a97f8cSStefano Zampini     PetscInt start_dof = queue_tip[i];
90734a97f8cSStefano Zampini     for (j=xadj[start_dof];j<xadj[start_dof+1];j++) {
90834a97f8cSStefano Zampini       PetscInt dof = adjncy[j];
90934a97f8cSStefano Zampini       if (!PetscBTLookup(touched,dof)) {
91034a97f8cSStefano Zampini         ierr = PetscBTSet(touched,dof);CHKERRQ(ierr);
91134a97f8cSStefano Zampini         queue_tip[n] = dof;
91234a97f8cSStefano Zampini         n++;
91334a97f8cSStefano Zampini       }
91434a97f8cSStefano Zampini     }
91534a97f8cSStefano Zampini   }
91634a97f8cSStefano Zampini   *n_added = n;
91734a97f8cSStefano Zampini   PetscFunctionReturn(0);
91834a97f8cSStefano Zampini }
919