xref: /petsc/src/ksp/pc/impls/bddc/bddcschurs.c (revision 2972d61b0bb859cde3c4726a4aafdfae2a88254b)
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;
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   }
535ec10c6aSStefano Zampini   if (!Bdense & !issym) {
543202ece2SStefano Zampini     ierr = MatDestroy(&Bd);CHKERRQ(ierr);
553202ece2SStefano Zampini   }
565ec10c6aSStefano Zampini 
575ec10c6aSStefano 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     }
635ec10c6aSStefano Zampini     ierr = MatMatMult(Cd, AinvBd, reuse, fill, S);CHKERRQ(ierr);
643202ece2SStefano Zampini     if (!Cdense) {
653202ece2SStefano Zampini       ierr = MatDestroy(&Cd);CHKERRQ(ierr);
663202ece2SStefano Zampini     }
675ec10c6aSStefano Zampini   } else {
685ec10c6aSStefano Zampini     ierr = MatTransposeMatMult(Bd, AinvBd, reuse, fill, S);CHKERRQ(ierr);
695ec10c6aSStefano Zampini     if (!Bdense) {
705ec10c6aSStefano Zampini       ierr = MatDestroy(&Bd);CHKERRQ(ierr);
715ec10c6aSStefano Zampini     }
725ec10c6aSStefano Zampini   }
735ec10c6aSStefano 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;
1065ec10c6aSStefano 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;
10906a4b1faSStefano Zampini   PetscScalar            *Bwork;
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 */
2055ec10c6aSStefano 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);
2105ec10c6aSStefano 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   }
264*2972d61bSStefano Zampini 
2653dc780c3SStefano Zampini   if (sub_schurs->n_subs) {
266883469d8SStefano Zampini     if (!sub_schurs->use_mumps) {
267883469d8SStefano Zampini       /* subsets in original and boundary numbering */
268883469d8SStefano Zampini       for (i=0;i<sub_schurs->n_subs;i++) {
269b96c3477SStefano Zampini         ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,sub_schurs->is_subs[i],&is_subset_B[i]);CHKERRQ(ierr);
270b1b3d7a2SStefano Zampini       }
271b1b3d7a2SStefano Zampini 
272b1b3d7a2SStefano Zampini       /* EE block */
273b1b3d7a2SStefano Zampini       for (i=0;i<sub_schurs->n_subs;i++) {
274b1b3d7a2SStefano Zampini         ierr = MatGetSubMatrix(A_BB,is_subset_B[i],is_subset_B[i],MAT_INITIAL_MATRIX,&AE_EE[i]);CHKERRQ(ierr);
275b1b3d7a2SStefano Zampini       }
276aa83b6aeSStefano Zampini 
277b1b3d7a2SStefano Zampini       /* IE block */
278b1b3d7a2SStefano Zampini       for (i=0;i<sub_schurs->n_subs;i++) {
279b1b3d7a2SStefano Zampini         ierr = MatGetSubMatrix(A_IB,is_I,is_subset_B[i],MAT_INITIAL_MATRIX,&AE_IE[i]);CHKERRQ(ierr);
280b1b3d7a2SStefano Zampini       }
281aa83b6aeSStefano Zampini 
282b1b3d7a2SStefano Zampini       /* EI block */
283b1b3d7a2SStefano Zampini       for (i=0;i<sub_schurs->n_subs;i++) {
284aa83b6aeSStefano Zampini         if (sub_schurs->is_hermitian) {
285aa83b6aeSStefano Zampini           ierr = MatCreateTranspose(AE_IE[i],&AE_EI[i]);CHKERRQ(ierr);
286aa83b6aeSStefano Zampini         } else {
287b1b3d7a2SStefano Zampini           ierr = MatGetSubMatrix(A_BI,is_subset_B[i],is_I,MAT_INITIAL_MATRIX,&AE_EI[i]);CHKERRQ(ierr);
288b1b3d7a2SStefano Zampini         }
289aa83b6aeSStefano Zampini       }
290b1b3d7a2SStefano Zampini 
291b1b3d7a2SStefano Zampini       /* setup Schur complements on subset */
292b1b3d7a2SStefano Zampini       for (i=0;i<sub_schurs->n_subs;i++) {
293b96c3477SStefano Zampini         ierr = MatDestroy(&sub_schurs->S_Ej[i]);CHKERRQ(ierr);
294b1b3d7a2SStefano Zampini         ierr = MatCreateSchurComplement(AE_II,AE_II,AE_IE[i],AE_EI[i],AE_EE[i],&sub_schurs->S_Ej[i]);CHKERRQ(ierr);
295a58a30b4SStefano Zampini         ierr = MatDestroy(&AE_EE[i]);CHKERRQ(ierr);
296a58a30b4SStefano Zampini         ierr = MatDestroy(&AE_IE[i]);CHKERRQ(ierr);
297a58a30b4SStefano Zampini         ierr = MatDestroy(&AE_EI[i]);CHKERRQ(ierr);
298b1b3d7a2SStefano Zampini         if (AE_II == A_II) { /* we can reuse the same ksp */
299b1b3d7a2SStefano Zampini           KSP ksp;
300b1b3d7a2SStefano Zampini           ierr = MatSchurComplementGetKSP(sub_schurs->S,&ksp);CHKERRQ(ierr);
301b1b3d7a2SStefano Zampini           ierr = MatSchurComplementSetKSP(sub_schurs->S_Ej[i],ksp);CHKERRQ(ierr);
302b1b3d7a2SStefano Zampini         } else { /* build new ksp object which inherits ksp and pc types from the original one */
303b1b3d7a2SStefano Zampini           KSP      origksp,schurksp;
304b1b3d7a2SStefano Zampini           PC       origpc,schurpc;
305b1b3d7a2SStefano Zampini           KSPType  ksp_type;
306b1b3d7a2SStefano Zampini           PCType   pc_type;
307b1b3d7a2SStefano Zampini           PetscInt n_internal;
308b1b3d7a2SStefano Zampini 
309b1b3d7a2SStefano Zampini           ierr = MatSchurComplementGetKSP(sub_schurs->S,&origksp);CHKERRQ(ierr);
310b1b3d7a2SStefano Zampini           ierr = MatSchurComplementGetKSP(sub_schurs->S_Ej[i],&schurksp);CHKERRQ(ierr);
311b1b3d7a2SStefano Zampini           ierr = KSPGetType(origksp,&ksp_type);CHKERRQ(ierr);
312b1b3d7a2SStefano Zampini           ierr = KSPSetType(schurksp,ksp_type);CHKERRQ(ierr);
313b1b3d7a2SStefano Zampini           ierr = KSPGetPC(schurksp,&schurpc);CHKERRQ(ierr);
314b1b3d7a2SStefano Zampini           ierr = KSPGetPC(origksp,&origpc);CHKERRQ(ierr);
315b1b3d7a2SStefano Zampini           ierr = PCGetType(origpc,&pc_type);CHKERRQ(ierr);
316b1b3d7a2SStefano Zampini           ierr = PCSetType(schurpc,pc_type);CHKERRQ(ierr);
317b1b3d7a2SStefano Zampini           ierr = ISGetSize(is_I,&n_internal);CHKERRQ(ierr);
318b1b3d7a2SStefano Zampini           if (n_internal) { /* UMFPACK gives error with 0 sized problems */
319b1b3d7a2SStefano Zampini             MatSolverPackage solver=NULL;
320b1b3d7a2SStefano Zampini             ierr = PCFactorGetMatSolverPackage(origpc,(const MatSolverPackage*)&solver);CHKERRQ(ierr);
321b1b3d7a2SStefano Zampini             if (solver) {
322b1b3d7a2SStefano Zampini               ierr = PCFactorSetMatSolverPackage(schurpc,solver);CHKERRQ(ierr);
323b1b3d7a2SStefano Zampini             }
324b1b3d7a2SStefano Zampini           }
325b1b3d7a2SStefano Zampini           ierr = KSPSetUp(schurksp);CHKERRQ(ierr);
326b1b3d7a2SStefano Zampini         }
327b1b3d7a2SStefano Zampini       }
328b1b3d7a2SStefano Zampini       /* free */
329b1b3d7a2SStefano Zampini       ierr = ISDestroy(&is_I);CHKERRQ(ierr);
330b1b3d7a2SStefano Zampini       ierr = MatDestroy(&AE_II);CHKERRQ(ierr);
331b1b3d7a2SStefano Zampini       for (i=0;i<sub_schurs->n_subs;i++) {
332b1b3d7a2SStefano Zampini         ierr = ISDestroy(&is_subset_B[i]);CHKERRQ(ierr);
333b1b3d7a2SStefano Zampini       }
334b1b3d7a2SStefano Zampini       ierr = PetscFree4(is_subset_B,AE_IE,AE_EI,AE_EE);CHKERRQ(ierr);
335883469d8SStefano Zampini #if defined(PETSC_HAVE_MUMPS)
336883469d8SStefano Zampini     } else {
337883469d8SStefano Zampini       Mat           A,F;
338883469d8SStefano Zampini       IS            is_A_all;
339883469d8SStefano Zampini       PetscInt      *idxs_schur,n_I;
340883469d8SStefano Zampini 
341883469d8SStefano Zampini       /* get working mat */
342883469d8SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_I_layer,&n_I);CHKERRQ(ierr);
343883469d8SStefano Zampini       ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size+n_I,all_local_idx_N,PETSC_COPY_VALUES,&is_A_all);CHKERRQ(ierr);
344b96c3477SStefano Zampini       ierr = MatGetSubMatrixUnsorted(sub_schurs->A,is_A_all,is_A_all,MAT_INITIAL_MATRIX,&A);CHKERRQ(ierr);
345883469d8SStefano Zampini       ierr = ISDestroy(&is_A_all);CHKERRQ(ierr);
346883469d8SStefano Zampini 
34708122e43SStefano Zampini       if (n_I) {
3489ab7bb16SStefano Zampini         if (sub_schurs->is_hermitian && sub_schurs->is_posdef) {
349883469d8SStefano Zampini           ierr = MatGetFactor(A,MATSOLVERMUMPS,MAT_FACTOR_CHOLESKY,&F);CHKERRQ(ierr);
350883469d8SStefano Zampini         } else {
351883469d8SStefano Zampini           ierr = MatGetFactor(A,MATSOLVERMUMPS,MAT_FACTOR_LU,&F);CHKERRQ(ierr);
352883469d8SStefano Zampini         }
353883469d8SStefano Zampini 
354883469d8SStefano Zampini         /* subsets ordered last */
355883469d8SStefano Zampini         ierr = PetscMalloc1(local_size,&idxs_schur);CHKERRQ(ierr);
356883469d8SStefano Zampini         for (i=0;i<local_size;i++) {
357883469d8SStefano Zampini           idxs_schur[i] = n_I+i+1;
358883469d8SStefano Zampini         }
359883469d8SStefano Zampini         ierr = MatMumpsSetSchurIndices(F,local_size,idxs_schur);CHKERRQ(ierr);
360883469d8SStefano Zampini         ierr = PetscFree(idxs_schur);CHKERRQ(ierr);
361883469d8SStefano Zampini 
362883469d8SStefano Zampini         /* factorization step */
3639ab7bb16SStefano Zampini         if (sub_schurs->is_hermitian && sub_schurs->is_posdef) {
364883469d8SStefano Zampini           ierr = MatCholeskyFactorSymbolic(F,A,NULL,NULL);CHKERRQ(ierr);
365883469d8SStefano Zampini           ierr = MatCholeskyFactorNumeric(F,A,NULL);CHKERRQ(ierr);
366883469d8SStefano Zampini         } else {
367883469d8SStefano Zampini           ierr = MatLUFactorSymbolic(F,A,NULL,NULL,NULL);CHKERRQ(ierr);
368883469d8SStefano Zampini           ierr = MatLUFactorNumeric(F,A,NULL);CHKERRQ(ierr);
369883469d8SStefano Zampini         }
370883469d8SStefano Zampini 
371883469d8SStefano Zampini         /* get explicit Schur Complement computed during numeric factorization */
372883469d8SStefano Zampini         ierr = MatMumpsGetSchurComplement(F,&S_all);CHKERRQ(ierr);
373883469d8SStefano Zampini         ierr = MatDestroy(&F);CHKERRQ(ierr);
37408122e43SStefano Zampini       } else {
37508122e43SStefano Zampini         ierr = MatConvert(A,MATSEQDENSE,MAT_INITIAL_MATRIX,&S_all);CHKERRQ(ierr);
37608122e43SStefano Zampini       }
37708122e43SStefano Zampini       ierr = MatDestroy(&A);CHKERRQ(ierr);
378883469d8SStefano Zampini #endif
379b1b3d7a2SStefano Zampini     }
3803dc780c3SStefano Zampini   } else {
3813dc780c3SStefano Zampini     ierr = PetscFree(nnz);CHKERRQ(ierr);
3823dc780c3SStefano Zampini     ierr = PetscFree(all_local_idx_N);CHKERRQ(ierr);
3833dc780c3SStefano Zampini   }
3845db18549SStefano Zampini 
3855db18549SStefano Zampini   if (!sub_schurs->n_subs_seq_g) {
3865db18549SStefano Zampini     PetscFunctionReturn(0);
3875db18549SStefano Zampini   }
3885db18549SStefano Zampini 
3895db18549SStefano Zampini   /* subset indices in local boundary numbering */
390b96c3477SStefano Zampini   if (!sub_schurs->is_Ej_all) {
391883469d8SStefano Zampini     ierr = ISGlobalToLocalMappingApply(sub_schurs->BtoNmap,IS_GTOLM_DROP,local_size,all_local_idx_N+extra,&subset_size,all_local_idx_B);CHKERRQ(ierr);
3925db18549SStefano Zampini     if (subset_size != local_size) {
3935db18549SStefano Zampini       SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in sub_schurs serial (BtoNmap)! %d != %d\n",subset_size,local_size);
3945db18549SStefano Zampini     }
3955db18549SStefano Zampini     ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size,all_local_idx_B,PETSC_OWN_POINTER,&sub_schurs->is_Ej_all);CHKERRQ(ierr);
396b96c3477SStefano Zampini   }
3975db18549SStefano Zampini 
3983202ece2SStefano Zampini   /* Local matrix of all local Schur on subsets transposed */
399b96c3477SStefano Zampini   if (!sub_schurs->S_Ej_all) {
4005db18549SStefano Zampini     ierr = MatCreate(PETSC_COMM_SELF,&sub_schurs->S_Ej_all);CHKERRQ(ierr);
4015db18549SStefano Zampini     ierr = MatSetSizes(sub_schurs->S_Ej_all,PETSC_DECIDE,PETSC_DECIDE,local_size,local_size);CHKERRQ(ierr);
4025db18549SStefano Zampini     ierr = MatSetType(sub_schurs->S_Ej_all,MATAIJ);CHKERRQ(ierr);
4035db18549SStefano Zampini     ierr = MatSeqAIJSetPreallocation(sub_schurs->S_Ej_all,0,nnz);CHKERRQ(ierr);
404b96c3477SStefano Zampini   } else {
405b96c3477SStefano Zampini     ierr = MatZeroEntries(sub_schurs->S_Ej_all);CHKERRQ(ierr);
406b96c3477SStefano Zampini   }
407a1337663SStefano Zampini 
408a1337663SStefano Zampini   S_Ej_tilda_all = 0;
40908122e43SStefano Zampini   S_Ej_inv_all = 0;
41006a4b1faSStefano Zampini   Bwork = NULL;
411a58a30b4SStefano Zampini   pivots = NULL;
412*2972d61bSStefano Zampini   if (sub_schurs->n_subs && deluxe && compute_Stilda && !sub_schurs->is_hermitian) { /* workspace needed only for GETRI */
41308122e43SStefano Zampini     PetscScalar lwork;
41408122e43SStefano Zampini 
41508122e43SStefano Zampini     B_lwork = -1;
4165ec10c6aSStefano Zampini     ierr = PetscBLASIntCast(max_subset_size_seq,&B_N);CHKERRQ(ierr);
41708122e43SStefano Zampini     ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
41806a4b1faSStefano Zampini     PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,Bwork,&B_N,pivots,&lwork,&B_lwork,&B_ierr));
41908122e43SStefano Zampini     ierr = PetscFPTrapPop();CHKERRQ(ierr);
42008122e43SStefano Zampini     if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to GETRI Lapack routine %d",(int)B_ierr);
42108122e43SStefano Zampini     ierr = PetscBLASIntCast((PetscInt)lwork,&B_lwork);CHKERRQ(ierr);
42206a4b1faSStefano Zampini     ierr = PetscMalloc2(B_lwork,&Bwork,max_subset_size_seq,&pivots);CHKERRQ(ierr);
42308122e43SStefano Zampini   }
42408122e43SStefano Zampini 
42508122e43SStefano Zampini   ierr = PetscBTMemzero(sub_schurs->n_subs,sub_schurs->computed_Stilda_subs);CHKERRQ(ierr);
426883469d8SStefano Zampini   if (!sub_schurs->use_mumps) {
4275ec10c6aSStefano Zampini     if (sub_schurs->n_subs_seq) {
4285ec10c6aSStefano Zampini       PetscScalar *work;
4295db18549SStefano Zampini       PetscInt    *dummy_idx;
4305db18549SStefano Zampini 
4315db18549SStefano Zampini       /* Work arrays */
4325ec10c6aSStefano Zampini       ierr = PetscMalloc2(max_subset_size_seq,&dummy_idx,max_subset_size_seq*max_subset_size_seq,&work);CHKERRQ(ierr);
4335ec10c6aSStefano Zampini 
4343202ece2SStefano Zampini       /* Loop on local problems end compute Schur complements explicitly */
4355db18549SStefano Zampini       local_size = 0;
4365db18549SStefano Zampini       for (i=0;i<sub_schurs->n_subs_seq;i++) {
4373202ece2SStefano Zampini         Mat       S_Ej_expl;
4383202ece2SStefano Zampini         PetscInt  j,lpi = sub_schurs->index_sequential[i];
4393202ece2SStefano Zampini         PetscBool Sdense;
4405db18549SStefano Zampini 
441b96c3477SStefano Zampini         ierr = ISGetLocalSize(sub_schurs->is_subs[lpi],&subset_size);CHKERRQ(ierr);
4425ec10c6aSStefano Zampini         ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&S_Ej_expl);CHKERRQ(ierr);
443aa83b6aeSStefano Zampini         ierr = PCBDDCComputeExplicitSchur(sub_schurs->S_Ej[lpi],sub_schurs->is_hermitian,MAT_REUSE_MATRIX,&S_Ej_expl);CHKERRQ(ierr);
4443202ece2SStefano Zampini         ierr = PetscObjectTypeCompare((PetscObject)S_Ej_expl,MATSEQDENSE,&Sdense);CHKERRQ(ierr);
4453202ece2SStefano Zampini         if (Sdense) {
4465db18549SStefano Zampini           for (j=0;j<subset_size;j++) {
4475db18549SStefano Zampini             dummy_idx[j]=local_size+j;
4485db18549SStefano Zampini           }
4495ec10c6aSStefano Zampini           ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
4503202ece2SStefano Zampini         } else {
4513202ece2SStefano Zampini           SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not yet implemented for sparse matrices");
4525db18549SStefano Zampini         }
4533202ece2SStefano Zampini         local_size += subset_size;
45465d8bf0aSStefano Zampini         ierr = MatDestroy(&sub_schurs->S_Ej[lpi]);CHKERRQ(ierr);
45565d8bf0aSStefano Zampini         sub_schurs->S_Ej[lpi] = S_Ej_expl;
4563202ece2SStefano Zampini       }
4575ec10c6aSStefano Zampini       ierr = PetscFree2(dummy_idx,work);CHKERRQ(ierr);
4585ec10c6aSStefano Zampini     }
459883469d8SStefano Zampini   } else {
460883469d8SStefano Zampini     PetscInt    *dummy_idx,n_all;
4615ec10c6aSStefano Zampini     PetscScalar *work;
462883469d8SStefano Zampini 
463a1337663SStefano Zampini     if (compute_Stilda) {
464a1337663SStefano Zampini       ierr = MatCreate(PETSC_COMM_SELF,&S_Ej_tilda_all);CHKERRQ(ierr);
465a1337663SStefano Zampini       ierr = MatSetSizes(S_Ej_tilda_all,PETSC_DECIDE,PETSC_DECIDE,local_size,local_size);CHKERRQ(ierr);
466a1337663SStefano Zampini       ierr = MatSetType(S_Ej_tilda_all,MATAIJ);CHKERRQ(ierr);
467a1337663SStefano Zampini       ierr = MatSeqAIJSetPreallocation(S_Ej_tilda_all,0,nnz);CHKERRQ(ierr);
46808122e43SStefano Zampini       if (deluxe) {
46908122e43SStefano Zampini         ierr = MatCreate(PETSC_COMM_SELF,&S_Ej_inv_all);CHKERRQ(ierr);
47008122e43SStefano Zampini         ierr = MatSetSizes(S_Ej_inv_all,PETSC_DECIDE,PETSC_DECIDE,local_size,local_size);CHKERRQ(ierr);
47108122e43SStefano Zampini         ierr = MatSetType(S_Ej_inv_all,MATAIJ);CHKERRQ(ierr);
47208122e43SStefano Zampini         ierr = MatSeqAIJSetPreallocation(S_Ej_inv_all,0,nnz);CHKERRQ(ierr);
473a1337663SStefano Zampini       }
47408122e43SStefano Zampini     }
47508122e43SStefano Zampini 
4768275c6a8SStefano Zampini     if (S_all) { /* multilevel */
47708122e43SStefano Zampini       ierr = MatGetSize(S_all,&n_all,NULL);CHKERRQ(ierr);
4788275c6a8SStefano Zampini     } else {
4798275c6a8SStefano Zampini       n_all = 0;
4808275c6a8SStefano Zampini     }
48108122e43SStefano Zampini     local_size = 0;
48206a4b1faSStefano Zampini 
483883469d8SStefano Zampini     /* Work arrays */
48406a4b1faSStefano Zampini     if (compute_Stilda) {
48506a4b1faSStefano Zampini       ierr = PetscMalloc2(max_subset_size_seq,&dummy_idx,n_all*n_all,&work);CHKERRQ(ierr);
48606a4b1faSStefano Zampini     } else {
48706a4b1faSStefano Zampini       ierr = PetscMalloc2(max_subset_size_seq,&dummy_idx,0,&work);CHKERRQ(ierr);
48806a4b1faSStefano Zampini     }
489883469d8SStefano Zampini 
49065d8bf0aSStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
49165d8bf0aSStefano Zampini       IS  is_E;
49265d8bf0aSStefano Zampini       PetscScalar *vals;
49365d8bf0aSStefano Zampini       PetscInt j;
49465d8bf0aSStefano Zampini 
495b96c3477SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
496a1337663SStefano Zampini       for (j=0;j<subset_size;j++) {
497a1337663SStefano Zampini         dummy_idx[j]=local_size+j;
498a1337663SStefano Zampini       }
49965d8bf0aSStefano Zampini       ierr = ISCreateStride(PetscObjectComm((PetscObject)sub_schurs->is_I),subset_size,local_size,1,&is_E);CHKERRQ(ierr);
500b96c3477SStefano Zampini       if (sub_schurs->S_Ej[i]) {
501b96c3477SStefano Zampini         ierr = MatGetSubMatrix(S_all,is_E,is_E,MAT_REUSE_MATRIX,&sub_schurs->S_Ej[i]);CHKERRQ(ierr);
502b96c3477SStefano Zampini       } else {
503a1337663SStefano Zampini         ierr = MatGetSubMatrix(S_all,is_E,is_E,MAT_INITIAL_MATRIX,&sub_schurs->S_Ej[i]);CHKERRQ(ierr);
504b96c3477SStefano Zampini       }
505a1337663SStefano Zampini       ierr = MatDenseGetArray(sub_schurs->S_Ej[i],&vals);CHKERRQ(ierr);
506a1337663SStefano Zampini       ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,vals,INSERT_VALUES);CHKERRQ(ierr);
507a1337663SStefano Zampini       ierr = MatDenseRestoreArray(sub_schurs->S_Ej[i],&vals);CHKERRQ(ierr);
508a1337663SStefano Zampini 
5096da91a30SStefano Zampini       if (compute_Stilda && ((PetscBTLookup(sub_schurs->is_edge,i) && use_edges) || use_faces)) {
51065d8bf0aSStefano Zampini         Mat Stilda;
51106a4b1faSStefano Zampini 
51206a4b1faSStefano Zampini         ierr = PetscBTSet(sub_schurs->computed_Stilda_subs,i);CHKERRQ(ierr);
51365d8bf0aSStefano Zampini         if (sub_schurs->n_subs == 1) {
51465d8bf0aSStefano Zampini           ierr = PetscObjectReference((PetscObject)sub_schurs->S_Ej[i]);CHKERRQ(ierr);
51565d8bf0aSStefano Zampini           Stilda = sub_schurs->S_Ej[i];
51665d8bf0aSStefano Zampini         } else {
51765d8bf0aSStefano Zampini           KSP         ksp;
51865d8bf0aSStefano Zampini           PC          pc;
51965d8bf0aSStefano Zampini           Mat         S_EF,S_FE,S_FF,Stilda_impl;
52065d8bf0aSStefano Zampini           IS          is_F;
5212b973097SStefano Zampini           PetscScalar eps=1.e-8;
52206a4b1faSStefano Zampini           PetscInt    nc,cum;
5232b973097SStefano Zampini           PetscBool   chop=PETSC_FALSE;
52465d8bf0aSStefano Zampini 
52565d8bf0aSStefano Zampini           ierr = ISComplement(is_E,0,n_all,&is_F);CHKERRQ(ierr);
52606a4b1faSStefano Zampini           nc   = n_all - subset_size;
52706a4b1faSStefano Zampini           cum  = 0;
52806a4b1faSStefano Zampini           ierr = MatCreateSeqDense(PETSC_COMM_SELF,nc,nc,work+cum,&S_FF);CHKERRQ(ierr);
529*2972d61bSStefano Zampini           cum += nc*nc;
53006a4b1faSStefano Zampini           ierr = MatCreateSeqDense(PETSC_COMM_SELF,nc,subset_size,work+cum,&S_FE);CHKERRQ(ierr);
531*2972d61bSStefano Zampini           cum += nc*subset_size;
53206a4b1faSStefano Zampini           ierr = MatGetSubMatrix(S_all,is_F,is_F,MAT_REUSE_MATRIX,&S_FF);CHKERRQ(ierr);
53306a4b1faSStefano Zampini           ierr = MatGetSubMatrix(S_all,is_F,is_E,MAT_REUSE_MATRIX,&S_FE);CHKERRQ(ierr);
5345ec10c6aSStefano Zampini           if (sub_schurs->is_hermitian) { /* just need a placeholder for S_EF */
53506a4b1faSStefano Zampini             ierr = MatCreateTranspose(S_FE,&S_EF);CHKERRQ(ierr);
5365ec10c6aSStefano Zampini           } else {
537*2972d61bSStefano Zampini             ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,nc,work+cum,&S_EF);CHKERRQ(ierr);
538*2972d61bSStefano Zampini             cum += nc*subset_size;
53906a4b1faSStefano Zampini             ierr = MatGetSubMatrix(S_all,is_E,is_F,MAT_REUSE_MATRIX,&S_EF);CHKERRQ(ierr);
5405ec10c6aSStefano Zampini           }
541a1337663SStefano Zampini           ierr = ISDestroy(&is_F);CHKERRQ(ierr);
5422b973097SStefano Zampini           if (chop) {
5432b973097SStefano Zampini             ierr = MatChop(S_FF,eps);CHKERRQ(ierr);
5442b973097SStefano Zampini             ierr = MatConvert(S_FF,MATAIJ,MAT_REUSE_MATRIX,&S_FF);CHKERRQ(ierr);
5452b973097SStefano Zampini           }
546a1337663SStefano Zampini           ierr = MatCreateSchurComplement(S_FF,S_FF,S_FE,S_EF,sub_schurs->S_Ej[i],&Stilda_impl);CHKERRQ(ierr);
547a1337663SStefano Zampini           ierr = MatDestroy(&S_FF);CHKERRQ(ierr);
548a1337663SStefano Zampini           ierr = MatDestroy(&S_FE);CHKERRQ(ierr);
549a1337663SStefano Zampini           ierr = MatDestroy(&S_EF);CHKERRQ(ierr);
55065d8bf0aSStefano Zampini           ierr = MatSchurComplementGetKSP(Stilda_impl,&ksp);CHKERRQ(ierr);
55165d8bf0aSStefano Zampini           ierr = KSPSetType(ksp,KSPPREONLY);CHKERRQ(ierr);
55265d8bf0aSStefano Zampini           ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr);
5539ab7bb16SStefano Zampini           if (sub_schurs->is_hermitian && sub_schurs->is_posdef) {
55465d8bf0aSStefano Zampini             ierr = PCSetType(pc,PCCHOLESKY);CHKERRQ(ierr);
55565d8bf0aSStefano Zampini           } else {
55665d8bf0aSStefano Zampini             ierr = PCSetType(pc,PCLU);CHKERRQ(ierr);
55765d8bf0aSStefano Zampini           }
5582b973097SStefano Zampini           if (!chop) {
55965d8bf0aSStefano Zampini             ierr = PCFactorSetUseInPlace(pc,PETSC_TRUE);CHKERRQ(ierr);
5602b973097SStefano Zampini           } else {
5612b973097SStefano Zampini             ierr = PCFactorSetMatSolverPackage(pc,MATSOLVERMUMPS);CHKERRQ(ierr);
5622b973097SStefano Zampini           }
56365d8bf0aSStefano Zampini           ierr = KSPSetUp(ksp);CHKERRQ(ierr);
56406a4b1faSStefano Zampini           ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work+cum,&Stilda);CHKERRQ(ierr);
5655ec10c6aSStefano Zampini           ierr = PCBDDCComputeExplicitSchur(Stilda_impl,sub_schurs->is_hermitian,MAT_REUSE_MATRIX,&Stilda);CHKERRQ(ierr);
56665d8bf0aSStefano Zampini           ierr = MatDestroy(&Stilda_impl);CHKERRQ(ierr);
56765d8bf0aSStefano Zampini         }
5682b973097SStefano Zampini /*
5692b973097SStefano Zampini         PetscViewerPushFormat(PETSC_VIEWER_STDOUT_SELF,PETSC_VIEWER_ASCII_MATLAB);
5702b973097SStefano Zampini         ierr = MatView(Stilda,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
5712b973097SStefano Zampini         PetscViewerPopFormat(PETSC_VIEWER_STDOUT_SELF);
5722b973097SStefano Zampini */
573b76a6013SStefano Zampini         ierr = MatDenseGetArray(Stilda,&vals);CHKERRQ(ierr);
57408122e43SStefano Zampini         if (deluxe) { /* when using deluxe scaling, we need (S_1^-1+S_2^-1)^-1 */
57508122e43SStefano Zampini           PetscScalar *vals2;
57608122e43SStefano Zampini 
57708122e43SStefano Zampini           ierr = MatDenseGetArray(sub_schurs->S_Ej[i],&vals2);CHKERRQ(ierr);
57808122e43SStefano Zampini           ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr);
57908122e43SStefano Zampini           ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
5809552c7c7SStefano Zampini           if (invert_Stildas) { /* Stildas can be singular */
581*2972d61bSStefano Zampini             if (!sub_schurs->is_hermitian) {
58208122e43SStefano Zampini               PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,vals,&B_N,pivots,&B_ierr));
58308122e43SStefano Zampini               if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
58406a4b1faSStefano Zampini               PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,vals,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
58508122e43SStefano Zampini               if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
586*2972d61bSStefano Zampini             } else {
587*2972d61bSStefano Zampini               PetscInt j,k;
588*2972d61bSStefano Zampini 
589*2972d61bSStefano Zampini               PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,vals,&B_N,&B_ierr));
590*2972d61bSStefano Zampini               if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
591*2972d61bSStefano Zampini               PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,vals,&B_N,&B_ierr));
592*2972d61bSStefano Zampini               if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
593*2972d61bSStefano Zampini               for (j=0;j<B_N;j++) {
594*2972d61bSStefano Zampini                 for (k=j+1;k<B_N;k++) {
595*2972d61bSStefano Zampini                   vals[k*B_N+j] = vals[j*B_N+k];
5969552c7c7SStefano Zampini                 }
597*2972d61bSStefano Zampini               }
598*2972d61bSStefano Zampini             }
599*2972d61bSStefano Zampini           }
600*2972d61bSStefano Zampini           if (!sub_schurs->is_hermitian) {
60108122e43SStefano Zampini             PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,vals2,&B_N,pivots,&B_ierr));
60208122e43SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
60306a4b1faSStefano Zampini             PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,vals2,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
60408122e43SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
605*2972d61bSStefano Zampini           } else {
606*2972d61bSStefano Zampini             PetscInt j,k;
607*2972d61bSStefano Zampini 
608*2972d61bSStefano Zampini             PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,vals2,&B_N,&B_ierr));
609*2972d61bSStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
610*2972d61bSStefano Zampini             PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,vals2,&B_N,&B_ierr));
611*2972d61bSStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
612*2972d61bSStefano Zampini             for (j=0;j<B_N;j++) {
613*2972d61bSStefano Zampini               for (k=j+1;k<B_N;k++) {
614*2972d61bSStefano Zampini                 vals2[k*B_N+j] = vals2[j*B_N+k];
615*2972d61bSStefano Zampini               }
616*2972d61bSStefano Zampini             }
617*2972d61bSStefano Zampini           }
61808122e43SStefano Zampini           ierr = PetscFPTrapPop();CHKERRQ(ierr);
61908122e43SStefano Zampini           ierr = MatSetValues(S_Ej_inv_all,subset_size,dummy_idx,subset_size,dummy_idx,vals2,INSERT_VALUES);CHKERRQ(ierr);
62008122e43SStefano Zampini           ierr = MatDenseRestoreArray(sub_schurs->S_Ej[i],&vals2);CHKERRQ(ierr);
62108122e43SStefano Zampini         }
622a1337663SStefano Zampini         ierr = MatSetValues(S_Ej_tilda_all,subset_size,dummy_idx,subset_size,dummy_idx,vals,INSERT_VALUES);CHKERRQ(ierr);
623b76a6013SStefano Zampini         ierr = MatDenseRestoreArray(Stilda,&vals);CHKERRQ(ierr);
624b76a6013SStefano Zampini         ierr = MatDestroy(&Stilda);CHKERRQ(ierr);
62565d8bf0aSStefano Zampini       }
62665d8bf0aSStefano Zampini       ierr = ISDestroy(&is_E);CHKERRQ(ierr);
627883469d8SStefano Zampini       local_size += subset_size;
628883469d8SStefano Zampini     }
6295ec10c6aSStefano Zampini     ierr = PetscFree2(dummy_idx,work);CHKERRQ(ierr);
6305db18549SStefano Zampini   }
631a1337663SStefano Zampini   ierr = PetscFree(nnz);CHKERRQ(ierr);
632a1337663SStefano Zampini   ierr = MatDestroy(&S_all);CHKERRQ(ierr);
6335db18549SStefano Zampini   ierr = MatAssemblyBegin(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
6345db18549SStefano Zampini   ierr = MatAssemblyEnd(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
6358275c6a8SStefano Zampini   if (S_Ej_tilda_all) {
636a1337663SStefano Zampini     ierr = MatAssemblyBegin(S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
637a1337663SStefano Zampini     ierr = MatAssemblyEnd(S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
6388275c6a8SStefano Zampini   }
6398275c6a8SStefano Zampini   if (S_Ej_inv_all) {
64008122e43SStefano Zampini     ierr = MatAssemblyBegin(S_Ej_inv_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
64108122e43SStefano Zampini     ierr = MatAssemblyEnd(S_Ej_inv_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
64208122e43SStefano Zampini   }
643a1337663SStefano Zampini 
6445db18549SStefano Zampini   /* Global matrix of all assembled Schur on subsets */
645883469d8SStefano 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);
6465db18549SStefano Zampini   ierr = ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)sub_schurs->l2gmap),1,local_size,all_local_idx_G,PETSC_COPY_VALUES,&l2gmap_subsets);CHKERRQ(ierr);
6475db18549SStefano Zampini   ierr = MatCreateIS(PetscObjectComm((PetscObject)sub_schurs->l2gmap),1,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size,l2gmap_subsets,&work_mat);CHKERRQ(ierr);
6485db18549SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&l2gmap_subsets);CHKERRQ(ierr);
6495db18549SStefano Zampini   ierr = MatISSetLocalMat(work_mat,sub_schurs->S_Ej_all);CHKERRQ(ierr);
6503927de2eSStefano Zampini   ierr = MatCreate(PetscObjectComm((PetscObject)work_mat),&global_schur_subsets);CHKERRQ(ierr);
6513927de2eSStefano Zampini   ierr = MatSetSizes(global_schur_subsets,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size);CHKERRQ(ierr);
6523927de2eSStefano Zampini   ierr = MatSetType(global_schur_subsets,MATMPIAIJ);CHKERRQ(ierr);
6533927de2eSStefano Zampini   ierr = MatISSetMPIXAIJPreallocation_Private(work_mat,global_schur_subsets,PETSC_TRUE);CHKERRQ(ierr);
6543927de2eSStefano Zampini   ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
6555db18549SStefano Zampini   /* Get local part of (\sum_j S_Ej) */
656883469d8SStefano Zampini   ierr = PetscFree(all_local_idx_N);CHKERRQ(ierr);
657883469d8SStefano Zampini   ierr = ISCreateGeneral(PetscObjectComm((PetscObject)sub_schurs->l2gmap),local_size,all_local_idx_G,PETSC_OWN_POINTER,&temp_is);CHKERRQ(ierr);
658b96c3477SStefano Zampini   if (sub_schurs->sum_S_Ej_all) {
659b96c3477SStefano Zampini     ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,MAT_REUSE_MATRIX,&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
660b96c3477SStefano Zampini   } else {
661b96c3477SStefano Zampini     ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,MAT_INITIAL_MATRIX,&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
662b96c3477SStefano Zampini   }
66308122e43SStefano Zampini 
6648275c6a8SStefano Zampini   if (S_Ej_tilda_all) {
665a1337663SStefano Zampini     ierr = MatISSetLocalMat(work_mat,S_Ej_tilda_all);CHKERRQ(ierr);
666a1337663SStefano Zampini     ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
667b96c3477SStefano Zampini     if (sub_schurs->sum_S_Ej_tilda_all) {
668b96c3477SStefano Zampini       ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,MAT_REUSE_MATRIX,&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
669b96c3477SStefano Zampini     } else {
670b96c3477SStefano Zampini       ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,MAT_INITIAL_MATRIX,&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
671b96c3477SStefano Zampini     }
6728275c6a8SStefano Zampini   }
6738275c6a8SStefano Zampini   if (S_Ej_inv_all) {
67408122e43SStefano Zampini     PetscInt    cum;
67508122e43SStefano Zampini     PetscScalar *array,*array2;
67608122e43SStefano Zampini     ierr = MatISSetLocalMat(work_mat,S_Ej_inv_all);CHKERRQ(ierr);
67708122e43SStefano Zampini     ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
67808122e43SStefano Zampini     if (sub_schurs->sum_S_Ej_inv_all) {
67908122e43SStefano Zampini       ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,MAT_REUSE_MATRIX,&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
68008122e43SStefano Zampini     } else {
68108122e43SStefano Zampini       ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,MAT_INITIAL_MATRIX,&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
68208122e43SStefano Zampini     }
68308122e43SStefano Zampini     /* invert blocks of sum_S_Ej_inv_all */
68408122e43SStefano Zampini     ierr = MatSeqAIJGetArray(sub_schurs->sum_S_Ej_inv_all,&array);CHKERRQ(ierr);
68508122e43SStefano Zampini     ierr = MatSeqAIJGetArray(sub_schurs->sum_S_Ej_tilda_all,&array2);CHKERRQ(ierr);
68608122e43SStefano Zampini     cum = 0;
68708122e43SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
68808122e43SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
68908122e43SStefano Zampini       if (PetscBTLookup(sub_schurs->computed_Stilda_subs,i)) {
69008122e43SStefano Zampini         ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr);
69108122e43SStefano Zampini         ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
692*2972d61bSStefano Zampini         if (!sub_schurs->is_hermitian) {
69308122e43SStefano Zampini           PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,array+cum,&B_N,pivots,&B_ierr));
69408122e43SStefano Zampini           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
69506a4b1faSStefano Zampini           PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,array+cum,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
69608122e43SStefano Zampini           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
697*2972d61bSStefano Zampini         } else {
698*2972d61bSStefano Zampini           PetscInt j,k;
699*2972d61bSStefano Zampini 
700*2972d61bSStefano Zampini           PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,array+cum,&B_N,&B_ierr));
701*2972d61bSStefano Zampini           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
702*2972d61bSStefano Zampini           PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,array+cum,&B_N,&B_ierr));
703*2972d61bSStefano Zampini           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
704*2972d61bSStefano Zampini           for (j=0;j<B_N;j++) {
705*2972d61bSStefano Zampini             for (k=j+1;k<B_N;k++) {
706*2972d61bSStefano Zampini               array[k*B_N+j+cum] = array[j*B_N+k+cum];
707*2972d61bSStefano Zampini             }
708*2972d61bSStefano Zampini           }
709*2972d61bSStefano Zampini         }
7109552c7c7SStefano Zampini         if (invert_Stildas) { /* Stildas can be singular */
711*2972d61bSStefano Zampini           if (!sub_schurs->is_hermitian) {
71208122e43SStefano Zampini             PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,array2+cum,&B_N,pivots,&B_ierr));
71308122e43SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
71406a4b1faSStefano Zampini             PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,array2+cum,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
71508122e43SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
716*2972d61bSStefano Zampini           } else {
717*2972d61bSStefano Zampini             PetscInt j,k;
718*2972d61bSStefano Zampini 
719*2972d61bSStefano Zampini             PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,array2+cum,&B_N,&B_ierr));
720*2972d61bSStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
721*2972d61bSStefano Zampini             PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,array2+cum,&B_N,&B_ierr));
722*2972d61bSStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
723*2972d61bSStefano Zampini             for (j=0;j<B_N;j++) {
724*2972d61bSStefano Zampini               for (k=j+1;k<B_N;k++) {
725*2972d61bSStefano Zampini                 array2[k*B_N+j+cum] = array2[j*B_N+k+cum];
726*2972d61bSStefano Zampini               }
727*2972d61bSStefano Zampini             }
728*2972d61bSStefano Zampini           }
7299552c7c7SStefano Zampini         }
73008122e43SStefano Zampini         ierr = PetscFPTrapPop();CHKERRQ(ierr);
73108122e43SStefano Zampini       }
73208122e43SStefano Zampini       cum += subset_size*subset_size;
73308122e43SStefano Zampini     }
73408122e43SStefano Zampini     ierr = MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_inv_all,&array);CHKERRQ(ierr);
73508122e43SStefano Zampini     ierr = MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_tilda_all,&array2);CHKERRQ(ierr);
73608122e43SStefano Zampini   }
7373202ece2SStefano Zampini 
73806a4b1faSStefano Zampini   ierr = PetscFree2(Bwork,pivots);CHKERRQ(ierr);
739a1337663SStefano Zampini   ierr = MatDestroy(&global_schur_subsets);CHKERRQ(ierr);
740a1337663SStefano Zampini   ierr = MatDestroy(&S_Ej_tilda_all);CHKERRQ(ierr);
74108122e43SStefano Zampini   ierr = MatDestroy(&S_Ej_inv_all);CHKERRQ(ierr);
7423202ece2SStefano Zampini   ierr = MatDestroy(&work_mat);CHKERRQ(ierr);
7435db18549SStefano Zampini   ierr = ISDestroy(&temp_is);CHKERRQ(ierr);
744b1b3d7a2SStefano Zampini   PetscFunctionReturn(0);
745b1b3d7a2SStefano Zampini }
746b1b3d7a2SStefano Zampini 
747b1b3d7a2SStefano Zampini #undef __FUNCT__
748b1b3d7a2SStefano Zampini #define __FUNCT__ "PCBDDCSubSchursInit"
7495db18549SStefano Zampini PetscErrorCode PCBDDCSubSchursInit(PCBDDCSubSchurs sub_schurs, Mat A, Mat S, IS is_I, IS is_B, PCBDDCGraph graph, ISLocalToGlobalMapping BtoNmap, PetscInt seqthreshold)
750b1b3d7a2SStefano Zampini {
7519bb4a8caSStefano Zampini   IS                  *faces,*edges,*all_cc,vertices;
752b1b3d7a2SStefano Zampini   PetscInt            *index_sequential,*index_parallel;
753b1b3d7a2SStefano Zampini   PetscInt            *auxlocal_sequential,*auxlocal_parallel;
754b1b3d7a2SStefano Zampini   PetscInt            *auxglobal_sequential,*auxglobal_parallel;
755b96c3477SStefano Zampini   PetscInt            *auxmapping;
756b1b3d7a2SStefano Zampini   PetscInt            i,max_subset_size;
757b1b3d7a2SStefano Zampini   PetscInt            n_sequential_problems,n_local_sequential_problems,n_parallel_problems,n_local_parallel_problems;
758b1b3d7a2SStefano Zampini   PetscInt            n_faces,n_edges,n_all_cc;
759b1b3d7a2SStefano Zampini   PetscBool           is_sorted;
760b1b3d7a2SStefano Zampini   PetscErrorCode  ierr;
761b1b3d7a2SStefano Zampini 
762b1b3d7a2SStefano Zampini   PetscFunctionBegin;
763b1b3d7a2SStefano Zampini   ierr = ISSorted(is_I,&is_sorted);CHKERRQ(ierr);
764b1b3d7a2SStefano Zampini   if (!is_sorted) {
765b1b3d7a2SStefano Zampini     SETERRQ(PetscObjectComm((PetscObject)is_I),PETSC_ERR_PLIB,"IS for I dofs should be shorted");
766b1b3d7a2SStefano Zampini   }
767b1b3d7a2SStefano Zampini   ierr = ISSorted(is_B,&is_sorted);CHKERRQ(ierr);
768b1b3d7a2SStefano Zampini   if (!is_sorted) {
769b1b3d7a2SStefano Zampini     SETERRQ(PetscObjectComm((PetscObject)is_B),PETSC_ERR_PLIB,"IS for B dofs should be shorted");
770b1b3d7a2SStefano Zampini   }
771b1b3d7a2SStefano Zampini 
772b1b3d7a2SStefano Zampini   /* reset any previous data */
773b1b3d7a2SStefano Zampini   ierr = PCBDDCSubSchursReset(sub_schurs);CHKERRQ(ierr);
774b1b3d7a2SStefano Zampini 
775b1b3d7a2SStefano Zampini   /* get index sets for faces and edges */
7769bb4a8caSStefano Zampini   ierr = PCBDDCGraphGetCandidatesIS(graph,&n_faces,&faces,&n_edges,&edges,&vertices);CHKERRQ(ierr);
777b1b3d7a2SStefano Zampini   n_all_cc = n_faces+n_edges;
77808122e43SStefano Zampini   ierr = PetscBTCreate(n_all_cc,&sub_schurs->is_edge);CHKERRQ(ierr);
77908122e43SStefano Zampini   ierr = PetscBTCreate(n_all_cc,&sub_schurs->computed_Stilda_subs);CHKERRQ(ierr);
780b1b3d7a2SStefano Zampini   ierr = PetscMalloc1(n_all_cc,&all_cc);CHKERRQ(ierr);
781b1b3d7a2SStefano Zampini   for (i=0;i<n_faces;i++) {
782b1b3d7a2SStefano Zampini     all_cc[i] = faces[i];
783b1b3d7a2SStefano Zampini   }
784b1b3d7a2SStefano Zampini   for (i=0;i<n_edges;i++) {
785b1b3d7a2SStefano Zampini     all_cc[n_faces+i] = edges[i];
78608122e43SStefano Zampini     ierr = PetscBTSet(sub_schurs->is_edge,n_faces+i);CHKERRQ(ierr);
787b1b3d7a2SStefano Zampini   }
788b1b3d7a2SStefano Zampini   ierr = PetscFree(faces);CHKERRQ(ierr);
789b1b3d7a2SStefano Zampini   ierr = PetscFree(edges);CHKERRQ(ierr);
790b1b3d7a2SStefano Zampini 
791b1b3d7a2SStefano Zampini   /* map interface's subsets */
792b1b3d7a2SStefano Zampini   max_subset_size = 0;
793b1b3d7a2SStefano Zampini   for (i=0;i<n_all_cc;i++) {
794b1b3d7a2SStefano Zampini     PetscInt subset_size;
795b1b3d7a2SStefano Zampini     ierr = ISGetLocalSize(all_cc[i],&subset_size);CHKERRQ(ierr);
796b1b3d7a2SStefano Zampini     max_subset_size = PetscMax(max_subset_size,subset_size);
797b1b3d7a2SStefano Zampini   }
798b1b3d7a2SStefano Zampini   ierr = PetscMalloc1(max_subset_size,&auxmapping);CHKERRQ(ierr);
799b1b3d7a2SStefano Zampini   ierr = PetscMalloc2(graph->ncc,&auxlocal_sequential,
800b1b3d7a2SStefano Zampini                       graph->ncc,&auxlocal_parallel);CHKERRQ(ierr);
801b1b3d7a2SStefano Zampini   ierr = PetscMalloc2(graph->ncc,&index_sequential,
802b1b3d7a2SStefano Zampini                       graph->ncc,&index_parallel);CHKERRQ(ierr);
803b1b3d7a2SStefano Zampini 
804883469d8SStefano Zampini   /* if threshold is negative uses all sequential problems (possibly using MUMPS) */
805883469d8SStefano Zampini   sub_schurs->use_mumps = PETSC_FALSE;
806883469d8SStefano Zampini   if (seqthreshold < 0) {
807883469d8SStefano Zampini     seqthreshold = max_subset_size;
808883469d8SStefano Zampini #if defined(PETSC_HAVE_MUMPS)
809883469d8SStefano Zampini     sub_schurs->use_mumps = !!A;
810883469d8SStefano Zampini #endif
811883469d8SStefano Zampini   }
812883469d8SStefano Zampini 
813b1b3d7a2SStefano Zampini 
814b1b3d7a2SStefano Zampini   /* determine which problem has to be solved in parallel or sequentially */
815b1b3d7a2SStefano Zampini   n_local_sequential_problems = 0;
816b1b3d7a2SStefano Zampini   n_local_parallel_problems = 0;
817b1b3d7a2SStefano Zampini   for (i=0;i<n_all_cc;i++) {
818b1b3d7a2SStefano Zampini     PetscInt       subset_size,j,min_loc = 0;
819b1b3d7a2SStefano Zampini     const PetscInt *idxs;
820b1b3d7a2SStefano Zampini 
821b1b3d7a2SStefano Zampini     ierr = ISGetLocalSize(all_cc[i],&subset_size);CHKERRQ(ierr);
822b1b3d7a2SStefano Zampini     ierr = ISGetIndices(all_cc[i],&idxs);CHKERRQ(ierr);
823b1b3d7a2SStefano Zampini     ierr = ISLocalToGlobalMappingApply(graph->l2gmap,subset_size,idxs,auxmapping);CHKERRQ(ierr);
824b1b3d7a2SStefano Zampini     for (j=1;j<subset_size;j++) {
825b1b3d7a2SStefano Zampini       if (auxmapping[j]<auxmapping[min_loc]) {
826b1b3d7a2SStefano Zampini         min_loc = j;
827b1b3d7a2SStefano Zampini       }
828b1b3d7a2SStefano Zampini     }
829b1b3d7a2SStefano Zampini     if (subset_size > seqthreshold) {
830b1b3d7a2SStefano Zampini       index_parallel[n_local_parallel_problems] = i;
831b1b3d7a2SStefano Zampini       auxlocal_parallel[n_local_parallel_problems] = idxs[min_loc];
832b1b3d7a2SStefano Zampini       n_local_parallel_problems++;
833b1b3d7a2SStefano Zampini     } else {
834b1b3d7a2SStefano Zampini       index_sequential[n_local_sequential_problems] = i;
835b1b3d7a2SStefano Zampini       auxlocal_sequential[n_local_sequential_problems] = idxs[min_loc];
836b1b3d7a2SStefano Zampini       n_local_sequential_problems++;
837b1b3d7a2SStefano Zampini     }
838b1b3d7a2SStefano Zampini     ierr = ISRestoreIndices(all_cc[i],&idxs);CHKERRQ(ierr);
839b1b3d7a2SStefano Zampini   }
840b1b3d7a2SStefano Zampini 
841b1b3d7a2SStefano Zampini   /* Number parallel problems */
842b1b3d7a2SStefano Zampini   auxglobal_parallel = 0;
843b1b3d7a2SStefano Zampini   ierr = PCBDDCSubsetNumbering(PetscObjectComm((PetscObject)graph->l2gmap),graph->l2gmap,n_local_parallel_problems,auxlocal_parallel,PETSC_NULL,&n_parallel_problems,&auxglobal_parallel);CHKERRQ(ierr);
844b1b3d7a2SStefano Zampini 
845b1b3d7a2SStefano Zampini   /* Number sequential problems */
846b1b3d7a2SStefano Zampini   auxglobal_sequential = 0;
847b1b3d7a2SStefano Zampini   ierr = PCBDDCSubsetNumbering(PetscObjectComm((PetscObject)graph->l2gmap),graph->l2gmap,n_local_sequential_problems,auxlocal_sequential,PETSC_NULL,&n_sequential_problems,&auxglobal_sequential);CHKERRQ(ierr);
848b1b3d7a2SStefano Zampini 
849b1b3d7a2SStefano Zampini   /* update info in sub_schurs */
850aa83b6aeSStefano Zampini   if (A) {
8519ab7bb16SStefano Zampini     PetscBool isseqaij;
8529ab7bb16SStefano Zampini 
8539ab7bb16SStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
8549ab7bb16SStefano Zampini     if (isseqaij) {
8551e9c79c2SStefano Zampini       ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
8561e9c79c2SStefano Zampini       sub_schurs->A = A;
8579ab7bb16SStefano Zampini     } else { /* SeqBAIJ matrices does not support symmetry checking, SeqSBAIJ does not support MatPermute */
8589ab7bb16SStefano Zampini       ierr = MatConvert(A,MATSEQAIJ,MAT_INITIAL_MATRIX,&sub_schurs->A);CHKERRQ(ierr);
8599ab7bb16SStefano Zampini     }
8601e9c79c2SStefano Zampini   }
861b1b3d7a2SStefano Zampini   ierr = PetscObjectReference((PetscObject)S);CHKERRQ(ierr);
862b1b3d7a2SStefano Zampini   sub_schurs->S = S;
863b1b3d7a2SStefano Zampini   ierr = PetscObjectReference((PetscObject)is_I);CHKERRQ(ierr);
864b1b3d7a2SStefano Zampini   sub_schurs->is_I = is_I;
865b1b3d7a2SStefano Zampini   ierr = PetscObjectReference((PetscObject)is_B);CHKERRQ(ierr);
866b1b3d7a2SStefano Zampini   sub_schurs->is_B = is_B;
8675db18549SStefano Zampini   ierr = PetscObjectReference((PetscObject)graph->l2gmap);CHKERRQ(ierr);
8685db18549SStefano Zampini   sub_schurs->l2gmap = graph->l2gmap;
8695db18549SStefano Zampini   ierr = PetscObjectReference((PetscObject)BtoNmap);CHKERRQ(ierr);
8705db18549SStefano Zampini   sub_schurs->BtoNmap = BtoNmap;
871b1b3d7a2SStefano Zampini   sub_schurs->n_subs_seq = n_local_sequential_problems;
872b1b3d7a2SStefano Zampini   sub_schurs->n_subs_par = n_local_parallel_problems;
873b1b3d7a2SStefano Zampini   sub_schurs->n_subs_seq_g = n_sequential_problems;
874b1b3d7a2SStefano Zampini   sub_schurs->n_subs_par_g = n_parallel_problems;
875b1b3d7a2SStefano Zampini   sub_schurs->n_subs = sub_schurs->n_subs_seq + sub_schurs->n_subs_par;
876b1b3d7a2SStefano Zampini   sub_schurs->is_subs = all_cc;
8779bb4a8caSStefano Zampini   if (!sub_schurs->use_mumps) { /* for adaptive selection */
878b96c3477SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
879b96c3477SStefano Zampini       ierr = ISSort(sub_schurs->is_subs[i]);CHKERRQ(ierr);
880b96c3477SStefano Zampini     }
8819bb4a8caSStefano Zampini   }
8829bb4a8caSStefano Zampini   sub_schurs->is_Ej_com = vertices;
883b1b3d7a2SStefano Zampini   sub_schurs->index_sequential = index_sequential;
884b1b3d7a2SStefano Zampini   sub_schurs->index_parallel = index_parallel;
885b1b3d7a2SStefano Zampini   sub_schurs->auxglobal_sequential = auxglobal_sequential;
886b1b3d7a2SStefano Zampini   sub_schurs->auxglobal_parallel = auxglobal_parallel;
887b1b3d7a2SStefano Zampini 
888b96c3477SStefano Zampini 
889b96c3477SStefano Zampini   /* allocate space for schur complements */
890b76a6013SStefano Zampini   ierr = PetscCalloc1(sub_schurs->n_subs,&sub_schurs->S_Ej);CHKERRQ(ierr);
891b96c3477SStefano Zampini   sub_schurs->S_Ej_all = NULL;
892b96c3477SStefano Zampini   sub_schurs->sum_S_Ej_all = NULL;
89308122e43SStefano Zampini   sub_schurs->sum_S_Ej_inv_all = NULL;
894b96c3477SStefano Zampini   sub_schurs->sum_S_Ej_tilda_all = NULL;
895b96c3477SStefano Zampini   sub_schurs->is_Ej_all = NULL;
896b96c3477SStefano Zampini 
897b1b3d7a2SStefano Zampini   /* free workspace */
898b1b3d7a2SStefano Zampini   ierr = PetscFree(auxmapping);CHKERRQ(ierr);
899b1b3d7a2SStefano Zampini   ierr = PetscFree2(auxlocal_sequential,auxlocal_parallel);CHKERRQ(ierr);
900b1b3d7a2SStefano Zampini 
901b1b3d7a2SStefano Zampini   PetscFunctionReturn(0);
902b1b3d7a2SStefano Zampini }
903b1b3d7a2SStefano Zampini 
904b1b3d7a2SStefano Zampini #undef __FUNCT__
90534a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCSubSchursCreate"
90634a97f8cSStefano Zampini PetscErrorCode PCBDDCSubSchursCreate(PCBDDCSubSchurs *sub_schurs)
90734a97f8cSStefano Zampini {
90834a97f8cSStefano Zampini   PCBDDCSubSchurs schurs_ctx;
90934a97f8cSStefano Zampini   PetscErrorCode  ierr;
91034a97f8cSStefano Zampini 
91134a97f8cSStefano Zampini   PetscFunctionBegin;
91234a97f8cSStefano Zampini   ierr = PetscNew(&schurs_ctx);CHKERRQ(ierr);
9135ff63025SStefano Zampini   schurs_ctx->n_subs = 0;
91434a97f8cSStefano Zampini   *sub_schurs = schurs_ctx;
91534a97f8cSStefano Zampini   PetscFunctionReturn(0);
91634a97f8cSStefano Zampini }
91734a97f8cSStefano Zampini 
91834a97f8cSStefano Zampini #undef __FUNCT__
91934a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCSubSchursDestroy"
92034a97f8cSStefano Zampini PetscErrorCode PCBDDCSubSchursDestroy(PCBDDCSubSchurs *sub_schurs)
92134a97f8cSStefano Zampini {
92234a97f8cSStefano Zampini   PetscErrorCode ierr;
92334a97f8cSStefano Zampini 
92434a97f8cSStefano Zampini   PetscFunctionBegin;
92534a97f8cSStefano Zampini   ierr = PCBDDCSubSchursReset(*sub_schurs);CHKERRQ(ierr);
92634a97f8cSStefano Zampini   ierr = PetscFree(*sub_schurs);CHKERRQ(ierr);
92734a97f8cSStefano Zampini   PetscFunctionReturn(0);
92834a97f8cSStefano Zampini }
92934a97f8cSStefano Zampini 
93034a97f8cSStefano Zampini #undef __FUNCT__
93134a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCSubSchursReset"
93234a97f8cSStefano Zampini PetscErrorCode PCBDDCSubSchursReset(PCBDDCSubSchurs sub_schurs)
93334a97f8cSStefano Zampini {
93434a97f8cSStefano Zampini   PetscInt       i;
93534a97f8cSStefano Zampini   PetscErrorCode ierr;
93634a97f8cSStefano Zampini 
93734a97f8cSStefano Zampini   PetscFunctionBegin;
9381e9c79c2SStefano Zampini   ierr = MatDestroy(&sub_schurs->A);CHKERRQ(ierr);
939b1b3d7a2SStefano Zampini   ierr = MatDestroy(&sub_schurs->S);CHKERRQ(ierr);
940b1b3d7a2SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_I);CHKERRQ(ierr);
941b1b3d7a2SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_B);CHKERRQ(ierr);
9425db18549SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->l2gmap);CHKERRQ(ierr);
9435db18549SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->BtoNmap);CHKERRQ(ierr);
94441c3ba1bSStefano Zampini   ierr = MatDestroy(&sub_schurs->S_Ej_all);CHKERRQ(ierr);
94541c3ba1bSStefano Zampini   ierr = MatDestroy(&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
94608122e43SStefano Zampini   ierr = MatDestroy(&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
947a1337663SStefano Zampini   ierr = MatDestroy(&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
9485db18549SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_Ej_all);CHKERRQ(ierr);
9499bb4a8caSStefano Zampini   ierr = ISDestroy(&sub_schurs->is_Ej_com);CHKERRQ(ierr);
95008122e43SStefano Zampini   ierr = PetscBTDestroy(&sub_schurs->is_edge);CHKERRQ(ierr);
95108122e43SStefano Zampini   ierr = PetscBTDestroy(&sub_schurs->computed_Stilda_subs);CHKERRQ(ierr);
95234a97f8cSStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
953b1b3d7a2SStefano Zampini     ierr = ISDestroy(&sub_schurs->is_subs[i]);CHKERRQ(ierr);
95434a97f8cSStefano Zampini     ierr = MatDestroy(&sub_schurs->S_Ej[i]);CHKERRQ(ierr);
95534a97f8cSStefano Zampini   }
95668270318SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_I_layer);CHKERRQ(ierr);
9575ff63025SStefano Zampini   if (sub_schurs->n_subs) {
958b1b3d7a2SStefano Zampini     ierr = PetscFree(sub_schurs->is_subs);CHKERRQ(ierr);
959b76a6013SStefano Zampini     ierr = PetscFree(sub_schurs->S_Ej);CHKERRQ(ierr);
9603dc780c3SStefano Zampini   }
961b1b3d7a2SStefano Zampini   ierr = PetscFree2(sub_schurs->index_sequential,sub_schurs->index_parallel);CHKERRQ(ierr);
962b1b3d7a2SStefano Zampini   ierr = PetscFree(sub_schurs->auxglobal_sequential);CHKERRQ(ierr);
963b1b3d7a2SStefano Zampini   ierr = PetscFree(sub_schurs->auxglobal_parallel);CHKERRQ(ierr);
96434a97f8cSStefano Zampini   sub_schurs->n_subs = 0;
96534a97f8cSStefano Zampini   PetscFunctionReturn(0);
96634a97f8cSStefano Zampini }
96734a97f8cSStefano Zampini 
96834a97f8cSStefano Zampini #undef __FUNCT__
96934a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCAdjGetNextLayer_Private"
9702a155e38SStefano Zampini PETSC_STATIC_INLINE PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt* queue_tip,PetscInt n_prev,PetscBT touched,PetscInt* xadj,PetscInt* adjncy,PetscInt* n_added)
97134a97f8cSStefano Zampini {
97234a97f8cSStefano Zampini   PetscInt       i,j,n;
97334a97f8cSStefano Zampini   PetscErrorCode ierr;
97434a97f8cSStefano Zampini 
97534a97f8cSStefano Zampini   PetscFunctionBegin;
97634a97f8cSStefano Zampini   n = 0;
97734a97f8cSStefano Zampini   for (i=-n_prev;i<0;i++) {
97834a97f8cSStefano Zampini     PetscInt start_dof = queue_tip[i];
97934a97f8cSStefano Zampini     for (j=xadj[start_dof];j<xadj[start_dof+1];j++) {
98034a97f8cSStefano Zampini       PetscInt dof = adjncy[j];
98134a97f8cSStefano Zampini       if (!PetscBTLookup(touched,dof)) {
98234a97f8cSStefano Zampini         ierr = PetscBTSet(touched,dof);CHKERRQ(ierr);
98334a97f8cSStefano Zampini         queue_tip[n] = dof;
98434a97f8cSStefano Zampini         n++;
98534a97f8cSStefano Zampini       }
98634a97f8cSStefano Zampini     }
98734a97f8cSStefano Zampini   }
98834a97f8cSStefano Zampini   *n_added = n;
98934a97f8cSStefano Zampini   PetscFunctionReturn(0);
99034a97f8cSStefano Zampini }
991