xref: /petsc/src/ksp/pc/impls/bddc/bddcschurs.c (revision d26273578b8ea1258ff979474489f51e08422fe2)
134a97f8cSStefano Zampini #include <../src/ksp/pc/impls/bddc/bddc.h>
234a97f8cSStefano Zampini #include <../src/ksp/pc/impls/bddc/bddcprivate.h>
308122e43SStefano Zampini #include <petscblaslapack.h>
434a97f8cSStefano Zampini 
53202ece2SStefano Zampini PETSC_STATIC_INLINE PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt*,PetscInt,PetscBT,PetscInt*,PetscInt*,PetscInt*);
65ec10c6aSStefano Zampini static PetscErrorCode PCBDDCComputeExplicitSchur(Mat,PetscBool,MatReuse,Mat*);
73202ece2SStefano Zampini 
83202ece2SStefano Zampini #undef __FUNCT__
93202ece2SStefano Zampini #define __FUNCT__ "PCBDDCComputeExplicitSchur"
105ec10c6aSStefano Zampini static PetscErrorCode PCBDDCComputeExplicitSchur(Mat M, PetscBool issym, MatReuse reuse, Mat *S)
113202ece2SStefano Zampini {
123202ece2SStefano Zampini   Mat            B, C, D, Bd, Cd, AinvBd;
133202ece2SStefano Zampini   KSP            ksp;
143202ece2SStefano Zampini   PC             pc;
153202ece2SStefano Zampini   PetscBool      isLU, isILU, isCHOL, Bdense, Cdense;
163202ece2SStefano Zampini   PetscReal      fill = 2.0;
17f11841e3SStefano Zampini   PetscInt       n_I;
183202ece2SStefano Zampini   PetscMPIInt    size;
193202ece2SStefano Zampini   PetscErrorCode ierr;
203202ece2SStefano Zampini 
213202ece2SStefano Zampini   PetscFunctionBegin;
223202ece2SStefano Zampini   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)M),&size);CHKERRQ(ierr);
233202ece2SStefano Zampini   if (size != 1) {
243202ece2SStefano Zampini     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not for parallel matrices");
253202ece2SStefano Zampini   }
26f11841e3SStefano Zampini   if (reuse == MAT_REUSE_MATRIX) {
27f11841e3SStefano Zampini     PetscBool Sdense;
28f11841e3SStefano Zampini 
29f11841e3SStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)*S, MATSEQDENSE, &Sdense);CHKERRQ(ierr);
30f11841e3SStefano Zampini     if (!Sdense) {
31f11841e3SStefano Zampini       SETERRQ(PetscObjectComm((PetscObject)M),PETSC_ERR_SUP,"S should dense");
32f11841e3SStefano Zampini     }
33f11841e3SStefano Zampini   }
343202ece2SStefano Zampini   ierr = MatSchurComplementGetSubMatrices(M, NULL, NULL, &B, &C, &D);CHKERRQ(ierr);
353202ece2SStefano Zampini   ierr = MatSchurComplementGetKSP(M, &ksp);CHKERRQ(ierr);
363202ece2SStefano Zampini   ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr);
373202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) pc, PCLU, &isLU);CHKERRQ(ierr);
383202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) pc, PCILU, &isILU);CHKERRQ(ierr);
393202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) pc, PCCHOLESKY, &isCHOL);CHKERRQ(ierr);
403202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) B, MATSEQDENSE, &Bdense);CHKERRQ(ierr);
413202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) C, MATSEQDENSE, &Cdense);CHKERRQ(ierr);
42f11841e3SStefano Zampini   ierr = MatGetSize(B,&n_I,NULL);CHKERRQ(ierr);
43f11841e3SStefano Zampini   if (n_I) {
443202ece2SStefano Zampini     if (!Bdense) {
453202ece2SStefano Zampini       ierr = MatConvert(B, MATSEQDENSE, MAT_INITIAL_MATRIX, &Bd);CHKERRQ(ierr);
463202ece2SStefano Zampini     } else {
473202ece2SStefano Zampini       Bd = B;
483202ece2SStefano Zampini     }
493202ece2SStefano Zampini 
503202ece2SStefano Zampini     if (isLU || isILU || isCHOL) {
513202ece2SStefano Zampini       Mat fact;
523202ece2SStefano Zampini       ierr = KSPSetUp(ksp);CHKERRQ(ierr);
533202ece2SStefano Zampini       ierr = PCFactorGetMatrix(pc, &fact);CHKERRQ(ierr);
543202ece2SStefano Zampini       ierr = MatDuplicate(Bd, MAT_DO_NOT_COPY_VALUES, &AinvBd);CHKERRQ(ierr);
553202ece2SStefano Zampini       ierr = MatMatSolve(fact, Bd, AinvBd);CHKERRQ(ierr);
563202ece2SStefano Zampini     } else {
5707b1e237SStefano Zampini       PetscBool ex = PETSC_TRUE;
5807b1e237SStefano Zampini 
5907b1e237SStefano Zampini       if (ex) {
603202ece2SStefano Zampini         Mat Ainvd;
613202ece2SStefano Zampini 
623202ece2SStefano Zampini         ierr = PCComputeExplicitOperator(pc, &Ainvd);CHKERRQ(ierr);
633202ece2SStefano Zampini         ierr = MatMatMult(Ainvd, Bd, MAT_INITIAL_MATRIX, fill, &AinvBd);CHKERRQ(ierr);
643202ece2SStefano Zampini         ierr = MatDestroy(&Ainvd);CHKERRQ(ierr);
6507b1e237SStefano Zampini       } else {
6607b1e237SStefano Zampini         Vec         sol,rhs;
6707b1e237SStefano Zampini         PetscScalar *arrayrhs,*arraysol;
6807b1e237SStefano Zampini         PetscInt    i,nrhs,n;
6907b1e237SStefano Zampini 
7007b1e237SStefano Zampini         ierr = MatDuplicate(Bd, MAT_DO_NOT_COPY_VALUES, &AinvBd);CHKERRQ(ierr);
7107b1e237SStefano Zampini         ierr = MatGetSize(Bd,&n,&nrhs);CHKERRQ(ierr);
7207b1e237SStefano Zampini         ierr = MatDenseGetArray(Bd,&arrayrhs);CHKERRQ(ierr);
7307b1e237SStefano Zampini         ierr = MatDenseGetArray(AinvBd,&arraysol);CHKERRQ(ierr);
7407b1e237SStefano Zampini         ierr = KSPGetSolution(ksp,&sol);CHKERRQ(ierr);
7507b1e237SStefano Zampini         ierr = KSPGetRhs(ksp,&rhs);CHKERRQ(ierr);
7607b1e237SStefano Zampini         for (i=0;i<nrhs;i++) {
7707b1e237SStefano Zampini           ierr = VecPlaceArray(rhs,arrayrhs+i*n);CHKERRQ(ierr);
7807b1e237SStefano Zampini           ierr = VecPlaceArray(sol,arraysol+i*n);CHKERRQ(ierr);
7907b1e237SStefano Zampini           ierr = KSPSolve(ksp,rhs,sol);CHKERRQ(ierr);
8007b1e237SStefano Zampini           ierr = VecResetArray(rhs);CHKERRQ(ierr);
8107b1e237SStefano Zampini           ierr = VecResetArray(sol);CHKERRQ(ierr);
8207b1e237SStefano Zampini         }
8307b1e237SStefano Zampini         ierr = MatDenseRestoreArray(Bd,&arrayrhs);CHKERRQ(ierr);
8407b1e237SStefano Zampini         ierr = MatDenseRestoreArray(AinvBd,&arrayrhs);CHKERRQ(ierr);
8507b1e237SStefano Zampini       }
863202ece2SStefano Zampini     }
875ec10c6aSStefano Zampini     if (!Bdense & !issym) {
883202ece2SStefano Zampini       ierr = MatDestroy(&Bd);CHKERRQ(ierr);
893202ece2SStefano Zampini     }
905ec10c6aSStefano Zampini 
915ec10c6aSStefano Zampini     if (!issym) {
923202ece2SStefano Zampini       if (!Cdense) {
933202ece2SStefano Zampini         ierr = MatConvert(C, MATSEQDENSE, MAT_INITIAL_MATRIX, &Cd);CHKERRQ(ierr);
943202ece2SStefano Zampini       } else {
953202ece2SStefano Zampini         Cd = C;
963202ece2SStefano Zampini       }
975ec10c6aSStefano Zampini       ierr = MatMatMult(Cd, AinvBd, reuse, fill, S);CHKERRQ(ierr);
983202ece2SStefano Zampini       if (!Cdense) {
993202ece2SStefano Zampini         ierr = MatDestroy(&Cd);CHKERRQ(ierr);
1003202ece2SStefano Zampini       }
1015ec10c6aSStefano Zampini     } else {
1025ec10c6aSStefano Zampini       ierr = MatTransposeMatMult(Bd, AinvBd, reuse, fill, S);CHKERRQ(ierr);
1035ec10c6aSStefano Zampini       if (!Bdense) {
1045ec10c6aSStefano Zampini         ierr = MatDestroy(&Bd);CHKERRQ(ierr);
1055ec10c6aSStefano Zampini       }
1065ec10c6aSStefano Zampini     }
1075ec10c6aSStefano Zampini     ierr = MatDestroy(&AinvBd);CHKERRQ(ierr);
108f11841e3SStefano Zampini   }
1093202ece2SStefano Zampini 
1103202ece2SStefano Zampini   if (D) {
1113202ece2SStefano Zampini     Mat       Dd;
1123202ece2SStefano Zampini     PetscBool Ddense;
1133202ece2SStefano Zampini 
1143202ece2SStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)D,MATSEQDENSE,&Ddense);CHKERRQ(ierr);
1153202ece2SStefano Zampini     if (!Ddense) {
1163202ece2SStefano Zampini       ierr = MatConvert(D, MATSEQDENSE, MAT_INITIAL_MATRIX, &Dd);CHKERRQ(ierr);
1173202ece2SStefano Zampini     } else {
1183202ece2SStefano Zampini       Dd = D;
1193202ece2SStefano Zampini     }
120f11841e3SStefano Zampini     if (n_I) {
1213202ece2SStefano Zampini       ierr = MatAYPX(*S,-1.0,Dd,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
122f11841e3SStefano Zampini     } else {
123f11841e3SStefano Zampini       if (reuse == MAT_INITIAL_MATRIX) {
124f11841e3SStefano Zampini         ierr = MatDuplicate(Dd,MAT_COPY_VALUES,S);CHKERRQ(ierr);
125f11841e3SStefano Zampini       } else {
126f11841e3SStefano Zampini         ierr = MatCopy(Dd,*S,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
127f11841e3SStefano Zampini       }
128f11841e3SStefano Zampini     }
1293202ece2SStefano Zampini     if (!Ddense) {
1303202ece2SStefano Zampini       ierr = MatDestroy(&Dd);CHKERRQ(ierr);
1313202ece2SStefano Zampini     }
1323202ece2SStefano Zampini   } else {
1333202ece2SStefano Zampini     ierr = MatScale(*S,-1.0);CHKERRQ(ierr);
1343202ece2SStefano Zampini   }
1353202ece2SStefano Zampini   PetscFunctionReturn(0);
1363202ece2SStefano Zampini }
13734a97f8cSStefano Zampini 
13834a97f8cSStefano Zampini #undef __FUNCT__
1391580ed26SStefano Zampini #define __FUNCT__ "PCBDDCSubSchursSetUp"
1409552c7c7SStefano 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)
141b1b3d7a2SStefano Zampini {
142b1b3d7a2SStefano Zampini   Mat                    A_II,A_IB,A_BI,A_BB;
143b1b3d7a2SStefano Zampini   Mat                    AE_II,*AE_IE,*AE_EI,*AE_EE;
144*d2627357SStefano Zampini   Mat                    S_all,S_all_inv;
145*d2627357SStefano Zampini   Mat                    global_schur_subsets,work_mat;
14608122e43SStefano Zampini   Mat                    S_Ej_tilda_all,S_Ej_inv_all;
1475db18549SStefano Zampini   ISLocalToGlobalMapping l2gmap_subsets;
1485db18549SStefano Zampini   IS                     is_I,*is_subset_B,temp_is;
149d648f858SStefano Zampini   PetscInt               *nnz,*all_local_idx_G,*all_local_idx_N;
1505ec10c6aSStefano Zampini   PetscInt               i,subset_size,max_subset_size_seq;
151883469d8SStefano Zampini   PetscInt               extra,local_size,global_size;
15208122e43SStefano Zampini   PetscBLASInt           B_N,B_ierr,B_lwork,*pivots;
15306a4b1faSStefano Zampini   PetscScalar            *Bwork;
154b1b3d7a2SStefano Zampini   PetscErrorCode         ierr;
155b1b3d7a2SStefano Zampini 
156b1b3d7a2SStefano Zampini   PetscFunctionBegin;
157b1b3d7a2SStefano Zampini   /* get Schur complement matrices */
158883469d8SStefano Zampini   if (!sub_schurs->use_mumps) {
159f11841e3SStefano Zampini     PetscBool isseqaij;
160f11841e3SStefano Zampini 
161b7eb3628SStefano Zampini     if (compute_Stilda) {
1623dc780c3SStefano Zampini       SETERRQ(PetscObjectComm((PetscObject)sub_schurs->l2gmap),PETSC_ERR_SUP,"Adaptive selection of constraints requires MUMPS");
163b7eb3628SStefano Zampini     }
164b1b3d7a2SStefano Zampini     ierr = MatSchurComplementGetSubMatrices(sub_schurs->S,&A_II,NULL,&A_IB,&A_BI,&A_BB);CHKERRQ(ierr);
165f11841e3SStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)A_BB,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
166f11841e3SStefano Zampini     if (!isseqaij) {
167f11841e3SStefano Zampini       ierr = MatConvert(A_BB,MATSEQAIJ,MAT_REUSE_MATRIX,&A_BB);CHKERRQ(ierr);
168f11841e3SStefano Zampini       ierr = MatConvert(A_IB,MATSEQAIJ,MAT_REUSE_MATRIX,&A_IB);CHKERRQ(ierr);
169f11841e3SStefano Zampini       ierr = MatConvert(A_BI,MATSEQAIJ,MAT_REUSE_MATRIX,&A_BI);CHKERRQ(ierr);
170f11841e3SStefano Zampini     }
171b1b3d7a2SStefano Zampini     ierr = PetscMalloc4(sub_schurs->n_subs,&is_subset_B,
172b1b3d7a2SStefano Zampini                         sub_schurs->n_subs,&AE_IE,
173b1b3d7a2SStefano Zampini                         sub_schurs->n_subs,&AE_EI,
174b1b3d7a2SStefano Zampini                         sub_schurs->n_subs,&AE_EE);CHKERRQ(ierr);
175a58a30b4SStefano Zampini   } else {
176a58a30b4SStefano Zampini     is_subset_B = NULL;
177a58a30b4SStefano Zampini     AE_IE = NULL;
178a58a30b4SStefano Zampini     AE_EI = NULL;
179a58a30b4SStefano Zampini     AE_EE = NULL;
180b1b3d7a2SStefano Zampini   }
181b1b3d7a2SStefano Zampini 
182b1b3d7a2SStefano Zampini   /* determine interior problems */
183b96c3477SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_I_layer);CHKERRQ(ierr);
1843dc780c3SStefano Zampini   ierr = ISGetLocalSize(sub_schurs->is_I,&i);CHKERRQ(ierr);
1853dc780c3SStefano Zampini   if (nlayers >= 0 && i) { /* Interior problems can be different from the original one */
186b1b3d7a2SStefano Zampini     PetscBT                touched;
187b1b3d7a2SStefano Zampini     const PetscInt*        idx_B;
188b1b3d7a2SStefano Zampini     PetscInt               n_I,n_B,n_local_dofs,n_prev_added,j,layer,*local_numbering;
189b1b3d7a2SStefano Zampini 
1903dc780c3SStefano Zampini     if (xadj == NULL || adjncy == NULL) {
1913dc780c3SStefano Zampini       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cannot request layering without adjacency");
1923dc780c3SStefano Zampini     }
193b1b3d7a2SStefano Zampini     /* get sizes */
194b1b3d7a2SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr);
195b1b3d7a2SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_B,&n_B);CHKERRQ(ierr);
196b1b3d7a2SStefano Zampini 
197b1b3d7a2SStefano Zampini     ierr = PetscMalloc1(n_I+n_B,&local_numbering);CHKERRQ(ierr);
198b1b3d7a2SStefano Zampini     ierr = PetscBTCreate(n_I+n_B,&touched);CHKERRQ(ierr);
199b1b3d7a2SStefano Zampini     ierr = PetscBTMemzero(n_I+n_B,touched);CHKERRQ(ierr);
200b1b3d7a2SStefano Zampini 
201b1b3d7a2SStefano Zampini     /* all boundary dofs must be skipped when adding layers */
202b1b3d7a2SStefano Zampini     ierr = ISGetIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr);
203b1b3d7a2SStefano Zampini     for (j=0;j<n_B;j++) {
204b1b3d7a2SStefano Zampini       ierr = PetscBTSet(touched,idx_B[j]);CHKERRQ(ierr);
205b1b3d7a2SStefano Zampini     }
206b1b3d7a2SStefano Zampini     ierr = PetscMemcpy(local_numbering,idx_B,n_B*sizeof(PetscInt));CHKERRQ(ierr);
207b1b3d7a2SStefano Zampini     ierr = ISRestoreIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr);
208b1b3d7a2SStefano Zampini 
209b1b3d7a2SStefano Zampini     /* add prescribed number of layers of dofs */
210b1b3d7a2SStefano Zampini     n_local_dofs = n_B;
211b1b3d7a2SStefano Zampini     n_prev_added = n_B;
212b1b3d7a2SStefano Zampini     for (layer=0;layer<nlayers;layer++) {
213b1b3d7a2SStefano Zampini       PetscInt n_added;
214b1b3d7a2SStefano Zampini       if (n_local_dofs == n_I+n_B) break;
215b1b3d7a2SStefano Zampini       if (n_local_dofs > n_I+n_B) {
216b1b3d7a2SStefano 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);
217b1b3d7a2SStefano Zampini       }
218b1b3d7a2SStefano Zampini       ierr = PCBDDCAdjGetNextLayer_Private(local_numbering+n_local_dofs,n_prev_added,touched,xadj,adjncy,&n_added);CHKERRQ(ierr);
219b1b3d7a2SStefano Zampini       n_prev_added = n_added;
220b1b3d7a2SStefano Zampini       n_local_dofs += n_added;
221b1b3d7a2SStefano Zampini       if (!n_added) break;
222b1b3d7a2SStefano Zampini     }
223b1b3d7a2SStefano Zampini     ierr = PetscBTDestroy(&touched);CHKERRQ(ierr);
224b1b3d7a2SStefano Zampini 
225883469d8SStefano Zampini     /* IS for I layer dofs in original numbering */
22668270318SStefano 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);
227b1b3d7a2SStefano Zampini     ierr = PetscFree(local_numbering);CHKERRQ(ierr);
22868270318SStefano Zampini     ierr = ISSort(sub_schurs->is_I_layer);CHKERRQ(ierr);
229883469d8SStefano Zampini     /* IS for I layer dofs in I numbering */
230883469d8SStefano Zampini     if (!sub_schurs->use_mumps) {
231b1b3d7a2SStefano Zampini       ISLocalToGlobalMapping ItoNmap;
232b1b3d7a2SStefano Zampini       ierr = ISLocalToGlobalMappingCreateIS(sub_schurs->is_I,&ItoNmap);CHKERRQ(ierr);
23368270318SStefano Zampini       ierr = ISGlobalToLocalMappingApplyIS(ItoNmap,IS_GTOLM_DROP,sub_schurs->is_I_layer,&is_I);CHKERRQ(ierr);
234b1b3d7a2SStefano Zampini       ierr = ISLocalToGlobalMappingDestroy(&ItoNmap);CHKERRQ(ierr);
235b1b3d7a2SStefano Zampini 
236b1b3d7a2SStefano Zampini       /* II block */
237b1b3d7a2SStefano Zampini       ierr = MatGetSubMatrix(A_II,is_I,is_I,MAT_INITIAL_MATRIX,&AE_II);CHKERRQ(ierr);
238b1b3d7a2SStefano Zampini     }
239b1b3d7a2SStefano Zampini   } else {
240b1b3d7a2SStefano Zampini     PetscInt n_I;
241b1b3d7a2SStefano Zampini 
242b1b3d7a2SStefano Zampini     /* IS for I dofs in original numbering */
243b1b3d7a2SStefano Zampini     ierr = PetscObjectReference((PetscObject)sub_schurs->is_I);CHKERRQ(ierr);
24468270318SStefano Zampini     sub_schurs->is_I_layer = sub_schurs->is_I;
245b1b3d7a2SStefano Zampini 
246b1b3d7a2SStefano Zampini     /* IS for I dofs in I numbering (strided 1) */
247883469d8SStefano Zampini     if (!sub_schurs->use_mumps) {
248b1b3d7a2SStefano Zampini       ierr = ISGetSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr);
249b1b3d7a2SStefano Zampini       ierr = ISCreateStride(PetscObjectComm((PetscObject)sub_schurs->is_I),n_I,0,1,&is_I);CHKERRQ(ierr);
250b1b3d7a2SStefano Zampini 
251b1b3d7a2SStefano Zampini       /* II block is the same */
252b1b3d7a2SStefano Zampini       ierr = PetscObjectReference((PetscObject)A_II);CHKERRQ(ierr);
253b1b3d7a2SStefano Zampini       AE_II = A_II;
254b1b3d7a2SStefano Zampini     }
255b1b3d7a2SStefano Zampini   }
256883469d8SStefano Zampini   /* Get info on subset sizes and sum of all subsets sizes */
2575ec10c6aSStefano Zampini   max_subset_size_seq = 0;
258883469d8SStefano Zampini   local_size = 0;
259883469d8SStefano Zampini   for (i=0;i<sub_schurs->n_subs_seq;i++) {
260883469d8SStefano Zampini     PetscInt j = sub_schurs->index_sequential[i];
261b96c3477SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_subs[j],&subset_size);CHKERRQ(ierr);
2625ec10c6aSStefano Zampini     max_subset_size_seq = PetscMax(subset_size,max_subset_size_seq);
263883469d8SStefano Zampini     local_size += subset_size;
264883469d8SStefano Zampini   }
265883469d8SStefano Zampini 
266883469d8SStefano Zampini   /* Work arrays for local indices */
267883469d8SStefano Zampini   extra = 0;
268883469d8SStefano Zampini   if (sub_schurs->use_mumps) {
269883469d8SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_I_layer,&extra);CHKERRQ(ierr);
270883469d8SStefano Zampini   }
271883469d8SStefano Zampini   ierr = PetscMalloc1(local_size+extra,&all_local_idx_N);CHKERRQ(ierr);
272883469d8SStefano Zampini   if (extra) {
273883469d8SStefano Zampini     const PetscInt *idxs;
274883469d8SStefano Zampini     ierr = ISGetIndices(sub_schurs->is_I_layer,&idxs);CHKERRQ(ierr);
275883469d8SStefano Zampini     ierr = PetscMemcpy(all_local_idx_N,idxs,extra*sizeof(PetscInt));CHKERRQ(ierr);
276883469d8SStefano Zampini     ierr = ISRestoreIndices(sub_schurs->is_I_layer,&idxs);CHKERRQ(ierr);
277883469d8SStefano Zampini   }
278883469d8SStefano Zampini   ierr = PetscMalloc1(local_size,&nnz);CHKERRQ(ierr);
279883469d8SStefano Zampini 
280883469d8SStefano Zampini   /* Get local indices in local numbering */
281883469d8SStefano Zampini   local_size = 0;
282883469d8SStefano Zampini   for (i=0;i<sub_schurs->n_subs_seq;i++) {
283883469d8SStefano Zampini     PetscInt j;
284883469d8SStefano Zampini     const    PetscInt *idxs;
285883469d8SStefano Zampini 
286883469d8SStefano Zampini     PetscInt local_problem_index = sub_schurs->index_sequential[i];
287b96c3477SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_subs[local_problem_index],&subset_size);CHKERRQ(ierr);
288b96c3477SStefano Zampini     ierr = ISGetIndices(sub_schurs->is_subs[local_problem_index],&idxs);CHKERRQ(ierr);
289883469d8SStefano Zampini     /* subset indices in local numbering */
290883469d8SStefano Zampini     ierr = PetscMemcpy(all_local_idx_N+local_size+extra,idxs,subset_size*sizeof(PetscInt));CHKERRQ(ierr);
291b96c3477SStefano Zampini     ierr = ISRestoreIndices(sub_schurs->is_subs[local_problem_index],&idxs);CHKERRQ(ierr);
292883469d8SStefano Zampini     for (j=0;j<subset_size;j++) nnz[local_size+j] = subset_size;
293883469d8SStefano Zampini     local_size += subset_size;
294883469d8SStefano Zampini   }
295883469d8SStefano Zampini 
296883469d8SStefano Zampini   S_all = NULL;
297*d2627357SStefano Zampini   S_all_inv = NULL;
298*d2627357SStefano Zampini   S_Ej_tilda_all = NULL;
299*d2627357SStefano Zampini   S_Ej_inv_all = NULL;
3009ab7bb16SStefano Zampini   sub_schurs->is_hermitian = PETSC_FALSE;
3019ab7bb16SStefano Zampini   sub_schurs->is_posdef = PETSC_FALSE;
3029552c7c7SStefano Zampini   if (sub_schurs->A) {
3039ab7bb16SStefano Zampini     ierr = MatIsHermitian(sub_schurs->A,0.0,&sub_schurs->is_hermitian);CHKERRQ(ierr);
3049ab7bb16SStefano Zampini     if (sub_schurs->is_hermitian) {
3059ab7bb16SStefano Zampini       PetscScalar val;
3069ab7bb16SStefano Zampini       Vec         vec1,vec2;
3079ab7bb16SStefano Zampini 
3089ab7bb16SStefano Zampini       ierr = MatCreateVecs(sub_schurs->A,&vec1,&vec2);CHKERRQ(ierr);
3099ab7bb16SStefano Zampini       ierr = VecSetRandom(vec1,NULL);
3109ab7bb16SStefano Zampini       ierr = VecCopy(vec1,vec2);CHKERRQ(ierr);
3119ab7bb16SStefano Zampini       ierr = MatMult(sub_schurs->A,vec2,vec1);CHKERRQ(ierr);
3129ab7bb16SStefano Zampini       ierr = VecDot(vec1,vec2,&val);CHKERRQ(ierr);
3139ab7bb16SStefano Zampini       if (PetscRealPart(val) > 0. && PetscImaginaryPart(val) == 0.) sub_schurs->is_posdef = PETSC_TRUE;
3149ab7bb16SStefano Zampini       ierr = VecDestroy(&vec1);CHKERRQ(ierr);
3159ab7bb16SStefano Zampini       ierr = VecDestroy(&vec2);CHKERRQ(ierr);
3169ab7bb16SStefano Zampini     }
3179552c7c7SStefano Zampini   }
318*d2627357SStefano Zampini   Bwork = NULL;
319*d2627357SStefano Zampini   pivots = NULL;
320*d2627357SStefano Zampini   if (sub_schurs->n_subs && deluxe && compute_Stilda && !sub_schurs->is_hermitian) { /* workspace needed only for GETRI */
321*d2627357SStefano Zampini     PetscScalar lwork;
322*d2627357SStefano Zampini 
323*d2627357SStefano Zampini     B_lwork = -1;
324*d2627357SStefano Zampini     ierr = PetscBLASIntCast(local_size,&B_N);CHKERRQ(ierr);
325*d2627357SStefano Zampini     ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
326*d2627357SStefano Zampini     PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,Bwork,&B_N,pivots,&lwork,&B_lwork,&B_ierr));
327*d2627357SStefano Zampini     ierr = PetscFPTrapPop();CHKERRQ(ierr);
328*d2627357SStefano Zampini     if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to GETRI Lapack routine %d",(int)B_ierr);
329*d2627357SStefano Zampini     ierr = PetscBLASIntCast((PetscInt)PetscRealPart(lwork),&B_lwork);CHKERRQ(ierr);
330*d2627357SStefano Zampini     ierr = PetscMalloc2(B_lwork,&Bwork,B_N,&pivots);CHKERRQ(ierr);
331*d2627357SStefano Zampini   }
332*d2627357SStefano Zampini 
333*d2627357SStefano Zampini   /* prepare parallel matrices for summing up properly schurs on subsets */
334*d2627357SStefano 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);
335*d2627357SStefano Zampini   ierr = ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)sub_schurs->l2gmap),1,local_size,all_local_idx_G,PETSC_COPY_VALUES,&l2gmap_subsets);CHKERRQ(ierr);
336*d2627357SStefano Zampini   ierr = MatCreateIS(PetscObjectComm((PetscObject)sub_schurs->l2gmap),1,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size,l2gmap_subsets,&work_mat);CHKERRQ(ierr);
337*d2627357SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&l2gmap_subsets);CHKERRQ(ierr);
338*d2627357SStefano Zampini   ierr = MatCreate(PetscObjectComm((PetscObject)work_mat),&global_schur_subsets);CHKERRQ(ierr);
339*d2627357SStefano Zampini   ierr = MatSetSizes(global_schur_subsets,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size);CHKERRQ(ierr);
340*d2627357SStefano Zampini   ierr = MatSetType(global_schur_subsets,MATMPIAIJ);CHKERRQ(ierr);
3412972d61bSStefano Zampini 
3423dc780c3SStefano Zampini   if (sub_schurs->n_subs) {
343883469d8SStefano Zampini     if (!sub_schurs->use_mumps) {
344883469d8SStefano Zampini       /* subsets in original and boundary numbering */
345883469d8SStefano Zampini       for (i=0;i<sub_schurs->n_subs;i++) {
346b96c3477SStefano Zampini         ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,sub_schurs->is_subs[i],&is_subset_B[i]);CHKERRQ(ierr);
347b1b3d7a2SStefano Zampini       }
348b1b3d7a2SStefano Zampini 
349b1b3d7a2SStefano Zampini       /* EE block */
350b1b3d7a2SStefano Zampini       for (i=0;i<sub_schurs->n_subs;i++) {
351b1b3d7a2SStefano Zampini         ierr = MatGetSubMatrix(A_BB,is_subset_B[i],is_subset_B[i],MAT_INITIAL_MATRIX,&AE_EE[i]);CHKERRQ(ierr);
352b1b3d7a2SStefano Zampini       }
353aa83b6aeSStefano Zampini 
354b1b3d7a2SStefano Zampini       /* IE block */
355b1b3d7a2SStefano Zampini       for (i=0;i<sub_schurs->n_subs;i++) {
356b1b3d7a2SStefano Zampini         ierr = MatGetSubMatrix(A_IB,is_I,is_subset_B[i],MAT_INITIAL_MATRIX,&AE_IE[i]);CHKERRQ(ierr);
357b1b3d7a2SStefano Zampini       }
358aa83b6aeSStefano Zampini 
359b1b3d7a2SStefano Zampini       /* EI block */
360b1b3d7a2SStefano Zampini       for (i=0;i<sub_schurs->n_subs;i++) {
361aa83b6aeSStefano Zampini         if (sub_schurs->is_hermitian) {
362aa83b6aeSStefano Zampini           ierr = MatCreateTranspose(AE_IE[i],&AE_EI[i]);CHKERRQ(ierr);
363aa83b6aeSStefano Zampini         } else {
364b1b3d7a2SStefano Zampini           ierr = MatGetSubMatrix(A_BI,is_subset_B[i],is_I,MAT_INITIAL_MATRIX,&AE_EI[i]);CHKERRQ(ierr);
365b1b3d7a2SStefano Zampini         }
366aa83b6aeSStefano Zampini       }
367b1b3d7a2SStefano Zampini 
368b1b3d7a2SStefano Zampini       /* setup Schur complements on subset */
369b1b3d7a2SStefano Zampini       for (i=0;i<sub_schurs->n_subs;i++) {
370b96c3477SStefano Zampini         ierr = MatDestroy(&sub_schurs->S_Ej[i]);CHKERRQ(ierr);
371b1b3d7a2SStefano Zampini         ierr = MatCreateSchurComplement(AE_II,AE_II,AE_IE[i],AE_EI[i],AE_EE[i],&sub_schurs->S_Ej[i]);CHKERRQ(ierr);
372a58a30b4SStefano Zampini         ierr = MatDestroy(&AE_EE[i]);CHKERRQ(ierr);
373a58a30b4SStefano Zampini         ierr = MatDestroy(&AE_IE[i]);CHKERRQ(ierr);
374a58a30b4SStefano Zampini         ierr = MatDestroy(&AE_EI[i]);CHKERRQ(ierr);
375b1b3d7a2SStefano Zampini         if (AE_II == A_II) { /* we can reuse the same ksp */
376b1b3d7a2SStefano Zampini           KSP ksp;
377b1b3d7a2SStefano Zampini           ierr = MatSchurComplementGetKSP(sub_schurs->S,&ksp);CHKERRQ(ierr);
378b1b3d7a2SStefano Zampini           ierr = MatSchurComplementSetKSP(sub_schurs->S_Ej[i],ksp);CHKERRQ(ierr);
379b1b3d7a2SStefano Zampini         } else { /* build new ksp object which inherits ksp and pc types from the original one */
380b1b3d7a2SStefano Zampini           KSP      origksp,schurksp;
381b1b3d7a2SStefano Zampini           PC       origpc,schurpc;
382b1b3d7a2SStefano Zampini           KSPType  ksp_type;
383b1b3d7a2SStefano Zampini           PCType   pc_type;
384b1b3d7a2SStefano Zampini           PetscInt n_internal;
385b1b3d7a2SStefano Zampini 
386b1b3d7a2SStefano Zampini           ierr = MatSchurComplementGetKSP(sub_schurs->S,&origksp);CHKERRQ(ierr);
387b1b3d7a2SStefano Zampini           ierr = MatSchurComplementGetKSP(sub_schurs->S_Ej[i],&schurksp);CHKERRQ(ierr);
388b1b3d7a2SStefano Zampini           ierr = KSPGetType(origksp,&ksp_type);CHKERRQ(ierr);
389b1b3d7a2SStefano Zampini           ierr = KSPSetType(schurksp,ksp_type);CHKERRQ(ierr);
390b1b3d7a2SStefano Zampini           ierr = KSPGetPC(schurksp,&schurpc);CHKERRQ(ierr);
391b1b3d7a2SStefano Zampini           ierr = KSPGetPC(origksp,&origpc);CHKERRQ(ierr);
392b1b3d7a2SStefano Zampini           ierr = PCGetType(origpc,&pc_type);CHKERRQ(ierr);
393b1b3d7a2SStefano Zampini           ierr = PCSetType(schurpc,pc_type);CHKERRQ(ierr);
394b1b3d7a2SStefano Zampini           ierr = ISGetSize(is_I,&n_internal);CHKERRQ(ierr);
395b1b3d7a2SStefano Zampini           if (n_internal) { /* UMFPACK gives error with 0 sized problems */
396b1b3d7a2SStefano Zampini             MatSolverPackage solver=NULL;
397b1b3d7a2SStefano Zampini             ierr = PCFactorGetMatSolverPackage(origpc,(const MatSolverPackage*)&solver);CHKERRQ(ierr);
398b1b3d7a2SStefano Zampini             if (solver) {
399b1b3d7a2SStefano Zampini               ierr = PCFactorSetMatSolverPackage(schurpc,solver);CHKERRQ(ierr);
400b1b3d7a2SStefano Zampini             }
401b1b3d7a2SStefano Zampini           }
402b1b3d7a2SStefano Zampini           ierr = KSPSetUp(schurksp);CHKERRQ(ierr);
403b1b3d7a2SStefano Zampini         }
404b1b3d7a2SStefano Zampini       }
405b1b3d7a2SStefano Zampini       /* free */
406b1b3d7a2SStefano Zampini       ierr = ISDestroy(&is_I);CHKERRQ(ierr);
407b1b3d7a2SStefano Zampini       ierr = MatDestroy(&AE_II);CHKERRQ(ierr);
408b1b3d7a2SStefano Zampini       for (i=0;i<sub_schurs->n_subs;i++) {
409b1b3d7a2SStefano Zampini         ierr = ISDestroy(&is_subset_B[i]);CHKERRQ(ierr);
410b1b3d7a2SStefano Zampini       }
411b1b3d7a2SStefano Zampini       ierr = PetscFree4(is_subset_B,AE_IE,AE_EI,AE_EE);CHKERRQ(ierr);
412883469d8SStefano Zampini #if defined(PETSC_HAVE_MUMPS)
413883469d8SStefano Zampini     } else {
414883469d8SStefano Zampini       Mat           A,F;
415883469d8SStefano Zampini       IS            is_A_all;
416883469d8SStefano Zampini       PetscInt      *idxs_schur,n_I;
417883469d8SStefano Zampini 
418883469d8SStefano Zampini       /* get working mat */
419883469d8SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_I_layer,&n_I);CHKERRQ(ierr);
420883469d8SStefano Zampini       ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size+n_I,all_local_idx_N,PETSC_COPY_VALUES,&is_A_all);CHKERRQ(ierr);
421d648f858SStefano Zampini       ierr = MatGetSubMatrixUnsorted(sub_schurs->A,is_A_all,is_A_all,&A);CHKERRQ(ierr);
422883469d8SStefano Zampini       ierr = ISDestroy(&is_A_all);CHKERRQ(ierr);
423883469d8SStefano Zampini 
42408122e43SStefano Zampini       if (n_I) {
4259ab7bb16SStefano Zampini         if (sub_schurs->is_hermitian && sub_schurs->is_posdef) {
426883469d8SStefano Zampini           ierr = MatGetFactor(A,MATSOLVERMUMPS,MAT_FACTOR_CHOLESKY,&F);CHKERRQ(ierr);
427883469d8SStefano Zampini         } else {
428883469d8SStefano Zampini           ierr = MatGetFactor(A,MATSOLVERMUMPS,MAT_FACTOR_LU,&F);CHKERRQ(ierr);
429883469d8SStefano Zampini         }
430883469d8SStefano Zampini 
431883469d8SStefano Zampini         /* subsets ordered last */
432883469d8SStefano Zampini         ierr = PetscMalloc1(local_size,&idxs_schur);CHKERRQ(ierr);
433883469d8SStefano Zampini         for (i=0;i<local_size;i++) {
434883469d8SStefano Zampini           idxs_schur[i] = n_I+i+1;
435883469d8SStefano Zampini         }
436883469d8SStefano Zampini         ierr = MatMumpsSetSchurIndices(F,local_size,idxs_schur);CHKERRQ(ierr);
437883469d8SStefano Zampini         ierr = PetscFree(idxs_schur);CHKERRQ(ierr);
438883469d8SStefano Zampini 
439883469d8SStefano Zampini         /* factorization step */
4409ab7bb16SStefano Zampini         if (sub_schurs->is_hermitian && sub_schurs->is_posdef) {
441883469d8SStefano Zampini           ierr = MatCholeskyFactorSymbolic(F,A,NULL,NULL);CHKERRQ(ierr);
442883469d8SStefano Zampini           ierr = MatCholeskyFactorNumeric(F,A,NULL);CHKERRQ(ierr);
443883469d8SStefano Zampini         } else {
444883469d8SStefano Zampini           ierr = MatLUFactorSymbolic(F,A,NULL,NULL,NULL);CHKERRQ(ierr);
445883469d8SStefano Zampini           ierr = MatLUFactorNumeric(F,A,NULL);CHKERRQ(ierr);
446883469d8SStefano Zampini         }
447883469d8SStefano Zampini 
448883469d8SStefano Zampini         /* get explicit Schur Complement computed during numeric factorization */
449883469d8SStefano Zampini         ierr = MatMumpsGetSchurComplement(F,&S_all);CHKERRQ(ierr);
450883469d8SStefano Zampini         ierr = MatDestroy(&F);CHKERRQ(ierr);
45108122e43SStefano Zampini       } else {
45208122e43SStefano Zampini         ierr = MatConvert(A,MATSEQDENSE,MAT_INITIAL_MATRIX,&S_all);CHKERRQ(ierr);
45308122e43SStefano Zampini       }
45408122e43SStefano Zampini       ierr = MatDestroy(&A);CHKERRQ(ierr);
455883469d8SStefano Zampini #endif
456b1b3d7a2SStefano Zampini     }
4573dc780c3SStefano Zampini   } else {
4583dc780c3SStefano Zampini     ierr = PetscFree(nnz);CHKERRQ(ierr);
4593dc780c3SStefano Zampini   }
4605db18549SStefano Zampini 
461*d2627357SStefano Zampini   if (!sub_schurs->n_subs_seq_g) { /* need to get rid of the old parallel deluxe */
462*d2627357SStefano Zampini     ierr = PetscFree(all_local_idx_G);CHKERRQ(ierr);
463*d2627357SStefano Zampini     ierr = PetscFree(all_local_idx_N);CHKERRQ(ierr);
464*d2627357SStefano Zampini     ierr = PetscFree2(Bwork,pivots);CHKERRQ(ierr);
465*d2627357SStefano Zampini     ierr = MatDestroy(&global_schur_subsets);CHKERRQ(ierr);
466*d2627357SStefano Zampini     ierr = MatDestroy(&work_mat);CHKERRQ(ierr);
4675db18549SStefano Zampini     PetscFunctionReturn(0);
4685db18549SStefano Zampini   }
4695db18549SStefano Zampini 
4705db18549SStefano Zampini   /* subset indices in local boundary numbering */
471b96c3477SStefano Zampini   if (!sub_schurs->is_Ej_all) {
472d648f858SStefano Zampini     PetscInt *all_local_idx_B;
473d648f858SStefano Zampini 
474d648f858SStefano Zampini     ierr = PetscMalloc1(local_size,&all_local_idx_B);CHKERRQ(ierr);
475883469d8SStefano Zampini     ierr = ISGlobalToLocalMappingApply(sub_schurs->BtoNmap,IS_GTOLM_DROP,local_size,all_local_idx_N+extra,&subset_size,all_local_idx_B);CHKERRQ(ierr);
4765db18549SStefano Zampini     if (subset_size != local_size) {
4775db18549SStefano Zampini       SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in sub_schurs serial (BtoNmap)! %d != %d\n",subset_size,local_size);
4785db18549SStefano Zampini     }
4795db18549SStefano Zampini     ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size,all_local_idx_B,PETSC_OWN_POINTER,&sub_schurs->is_Ej_all);CHKERRQ(ierr);
480b96c3477SStefano Zampini   }
481*d2627357SStefano Zampini   ierr = PetscFree(all_local_idx_N);CHKERRQ(ierr);
4825db18549SStefano Zampini 
483*d2627357SStefano Zampini   /* Local matrix of all local Schur on subsets (transposed) */
484b96c3477SStefano Zampini   if (!sub_schurs->S_Ej_all) {
4855db18549SStefano Zampini     ierr = MatCreate(PETSC_COMM_SELF,&sub_schurs->S_Ej_all);CHKERRQ(ierr);
4865db18549SStefano Zampini     ierr = MatSetSizes(sub_schurs->S_Ej_all,PETSC_DECIDE,PETSC_DECIDE,local_size,local_size);CHKERRQ(ierr);
4875db18549SStefano Zampini     ierr = MatSetType(sub_schurs->S_Ej_all,MATAIJ);CHKERRQ(ierr);
4885db18549SStefano Zampini     ierr = MatSeqAIJSetPreallocation(sub_schurs->S_Ej_all,0,nnz);CHKERRQ(ierr);
489b96c3477SStefano Zampini   } else {
490b96c3477SStefano Zampini     ierr = MatZeroEntries(sub_schurs->S_Ej_all);CHKERRQ(ierr);
491b96c3477SStefano Zampini   }
492a1337663SStefano Zampini 
49308122e43SStefano Zampini   ierr = PetscBTMemzero(sub_schurs->n_subs,sub_schurs->computed_Stilda_subs);CHKERRQ(ierr);
494883469d8SStefano Zampini   if (!sub_schurs->use_mumps) {
4955ec10c6aSStefano Zampini     if (sub_schurs->n_subs_seq) {
4965ec10c6aSStefano Zampini       PetscScalar *work;
4975db18549SStefano Zampini       PetscInt    *dummy_idx;
4985db18549SStefano Zampini 
4995db18549SStefano Zampini       /* Work arrays */
5005ec10c6aSStefano Zampini       ierr = PetscMalloc2(max_subset_size_seq,&dummy_idx,max_subset_size_seq*max_subset_size_seq,&work);CHKERRQ(ierr);
5015ec10c6aSStefano Zampini 
5023202ece2SStefano Zampini       /* Loop on local problems end compute Schur complements explicitly */
5035db18549SStefano Zampini       local_size = 0;
5045db18549SStefano Zampini       for (i=0;i<sub_schurs->n_subs_seq;i++) {
5053202ece2SStefano Zampini         Mat       S_Ej_expl;
5063202ece2SStefano Zampini         PetscInt  j,lpi = sub_schurs->index_sequential[i];
5073202ece2SStefano Zampini         PetscBool Sdense;
5085db18549SStefano Zampini 
509b96c3477SStefano Zampini         ierr = ISGetLocalSize(sub_schurs->is_subs[lpi],&subset_size);CHKERRQ(ierr);
5105ec10c6aSStefano Zampini         ierr = MatCreateSeqDense(PETSC_COMM_SELF,subset_size,subset_size,work,&S_Ej_expl);CHKERRQ(ierr);
511aa83b6aeSStefano Zampini         ierr = PCBDDCComputeExplicitSchur(sub_schurs->S_Ej[lpi],sub_schurs->is_hermitian,MAT_REUSE_MATRIX,&S_Ej_expl);CHKERRQ(ierr);
5123202ece2SStefano Zampini         ierr = PetscObjectTypeCompare((PetscObject)S_Ej_expl,MATSEQDENSE,&Sdense);CHKERRQ(ierr);
5133202ece2SStefano Zampini         if (Sdense) {
5145db18549SStefano Zampini           for (j=0;j<subset_size;j++) {
5155db18549SStefano Zampini             dummy_idx[j]=local_size+j;
5165db18549SStefano Zampini           }
5175ec10c6aSStefano Zampini           ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,work,INSERT_VALUES);CHKERRQ(ierr);
5183202ece2SStefano Zampini         } else {
5193202ece2SStefano Zampini           SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not yet implemented for sparse matrices");
5205db18549SStefano Zampini         }
5213202ece2SStefano Zampini         local_size += subset_size;
52265d8bf0aSStefano Zampini         ierr = MatDestroy(&sub_schurs->S_Ej[lpi]);CHKERRQ(ierr);
52365d8bf0aSStefano Zampini         sub_schurs->S_Ej[lpi] = S_Ej_expl;
5243202ece2SStefano Zampini       }
5255ec10c6aSStefano Zampini       ierr = PetscFree2(dummy_idx,work);CHKERRQ(ierr);
5265ec10c6aSStefano Zampini     }
527883469d8SStefano Zampini   } else {
528883469d8SStefano Zampini     PetscInt    *dummy_idx,n_all;
5295ec10c6aSStefano Zampini     PetscScalar *work;
530883469d8SStefano Zampini 
531a1337663SStefano Zampini     if (compute_Stilda) {
532a1337663SStefano Zampini       ierr = MatCreate(PETSC_COMM_SELF,&S_Ej_tilda_all);CHKERRQ(ierr);
533a1337663SStefano Zampini       ierr = MatSetSizes(S_Ej_tilda_all,PETSC_DECIDE,PETSC_DECIDE,local_size,local_size);CHKERRQ(ierr);
534a1337663SStefano Zampini       ierr = MatSetType(S_Ej_tilda_all,MATAIJ);CHKERRQ(ierr);
535a1337663SStefano Zampini       ierr = MatSeqAIJSetPreallocation(S_Ej_tilda_all,0,nnz);CHKERRQ(ierr);
53608122e43SStefano Zampini       if (deluxe) {
53708122e43SStefano Zampini         ierr = MatCreate(PETSC_COMM_SELF,&S_Ej_inv_all);CHKERRQ(ierr);
53808122e43SStefano Zampini         ierr = MatSetSizes(S_Ej_inv_all,PETSC_DECIDE,PETSC_DECIDE,local_size,local_size);CHKERRQ(ierr);
53908122e43SStefano Zampini         ierr = MatSetType(S_Ej_inv_all,MATAIJ);CHKERRQ(ierr);
54008122e43SStefano Zampini         ierr = MatSeqAIJSetPreallocation(S_Ej_inv_all,0,nnz);CHKERRQ(ierr);
541a1337663SStefano Zampini       }
54208122e43SStefano Zampini     }
54308122e43SStefano Zampini 
5448275c6a8SStefano Zampini     if (S_all) { /* multilevel */
54508122e43SStefano Zampini       ierr = MatGetSize(S_all,&n_all,NULL);CHKERRQ(ierr);
5468275c6a8SStefano Zampini     } else {
5478275c6a8SStefano Zampini       n_all = 0;
5488275c6a8SStefano Zampini     }
54908122e43SStefano Zampini     local_size = 0;
55006a4b1faSStefano Zampini 
551883469d8SStefano Zampini     /* Work arrays */
55206a4b1faSStefano Zampini     if (compute_Stilda) {
55306a4b1faSStefano Zampini       ierr = PetscMalloc2(max_subset_size_seq,&dummy_idx,n_all*n_all,&work);CHKERRQ(ierr);
55406a4b1faSStefano Zampini     } else {
55506a4b1faSStefano Zampini       ierr = PetscMalloc2(max_subset_size_seq,&dummy_idx,0,&work);CHKERRQ(ierr);
55606a4b1faSStefano Zampini     }
557883469d8SStefano Zampini 
558*d2627357SStefano Zampini     if (compute_Stilda) {
559*d2627357SStefano Zampini       PetscScalar *vals;
560*d2627357SStefano Zampini       ierr = PetscBLASIntCast(n_all,&B_N);CHKERRQ(ierr);
561*d2627357SStefano Zampini       ierr = MatDuplicate(S_all,MAT_COPY_VALUES,&S_all_inv);CHKERRQ(ierr);
562*d2627357SStefano Zampini       ierr = MatDenseGetArray(S_all_inv,&vals);CHKERRQ(ierr);
563*d2627357SStefano Zampini       if (!sub_schurs->is_hermitian) {
564*d2627357SStefano Zampini         PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,vals,&B_N,pivots,&B_ierr));
565*d2627357SStefano Zampini         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
566*d2627357SStefano Zampini         PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,vals,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
567*d2627357SStefano Zampini         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
568*d2627357SStefano Zampini       } else {
569*d2627357SStefano Zampini         PetscInt j,k;
570*d2627357SStefano Zampini 
571*d2627357SStefano Zampini         PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,vals,&B_N,&B_ierr));
572*d2627357SStefano Zampini         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
573*d2627357SStefano Zampini         PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,vals,&B_N,&B_ierr));
574*d2627357SStefano Zampini         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
575*d2627357SStefano Zampini         for (j=0;j<B_N;j++) {
576*d2627357SStefano Zampini           for (k=j+1;k<B_N;k++) {
577*d2627357SStefano Zampini             vals[k*B_N+j] = vals[j*B_N+k];
578*d2627357SStefano Zampini           }
579*d2627357SStefano Zampini         }
580*d2627357SStefano Zampini       }
581*d2627357SStefano Zampini       ierr = MatDenseRestoreArray(S_all_inv,&vals);CHKERRQ(ierr);
582*d2627357SStefano Zampini     }
583*d2627357SStefano Zampini 
58465d8bf0aSStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
58565d8bf0aSStefano Zampini       IS  is_E;
58665d8bf0aSStefano Zampini       PetscScalar *vals;
58765d8bf0aSStefano Zampini       PetscInt j;
58865d8bf0aSStefano Zampini 
589b96c3477SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
590a1337663SStefano Zampini       for (j=0;j<subset_size;j++) {
591a1337663SStefano Zampini         dummy_idx[j]=local_size+j;
592a1337663SStefano Zampini       }
59365d8bf0aSStefano Zampini       ierr = ISCreateStride(PetscObjectComm((PetscObject)sub_schurs->is_I),subset_size,local_size,1,&is_E);CHKERRQ(ierr);
594b96c3477SStefano Zampini       if (sub_schurs->S_Ej[i]) {
595b96c3477SStefano Zampini         ierr = MatGetSubMatrix(S_all,is_E,is_E,MAT_REUSE_MATRIX,&sub_schurs->S_Ej[i]);CHKERRQ(ierr);
596b96c3477SStefano Zampini       } else {
597a1337663SStefano Zampini         ierr = MatGetSubMatrix(S_all,is_E,is_E,MAT_INITIAL_MATRIX,&sub_schurs->S_Ej[i]);CHKERRQ(ierr);
598b96c3477SStefano Zampini       }
599a1337663SStefano Zampini       ierr = MatDenseGetArray(sub_schurs->S_Ej[i],&vals);CHKERRQ(ierr);
600a1337663SStefano Zampini       ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,vals,INSERT_VALUES);CHKERRQ(ierr);
601a1337663SStefano Zampini       ierr = MatDenseRestoreArray(sub_schurs->S_Ej[i],&vals);CHKERRQ(ierr);
602a1337663SStefano Zampini 
603*d2627357SStefano Zampini       if (compute_Stilda && ((PetscBTLookup(sub_schurs->is_edge,i) && use_edges) || (!PetscBTLookup(sub_schurs->is_edge,i) && use_faces))) {
60465d8bf0aSStefano Zampini         Mat Stilda;
60506a4b1faSStefano Zampini 
60606a4b1faSStefano Zampini         ierr = PetscBTSet(sub_schurs->computed_Stilda_subs,i);CHKERRQ(ierr);
607*d2627357SStefano Zampini         ierr = MatGetSubMatrix(S_all_inv,is_E,is_E,MAT_INITIAL_MATRIX,&Stilda);CHKERRQ(ierr);
608b76a6013SStefano Zampini         ierr = MatDenseGetArray(Stilda,&vals);CHKERRQ(ierr);
60908122e43SStefano Zampini         if (deluxe) { /* when using deluxe scaling, we need (S_1^-1+S_2^-1)^-1 */
61008122e43SStefano Zampini           PetscScalar *vals2;
61108122e43SStefano Zampini 
61208122e43SStefano Zampini           ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr);
61308122e43SStefano Zampini           ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
614*d2627357SStefano Zampini           ierr = MatDenseGetArray(sub_schurs->S_Ej[i],&vals2);CHKERRQ(ierr);
6152972d61bSStefano Zampini           if (!sub_schurs->is_hermitian) {
61608122e43SStefano Zampini             PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,vals2,&B_N,pivots,&B_ierr));
61708122e43SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
61806a4b1faSStefano Zampini             PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,vals2,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
61908122e43SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
6202972d61bSStefano Zampini           } else {
6212972d61bSStefano Zampini             PetscInt j,k;
6222972d61bSStefano Zampini 
6232972d61bSStefano Zampini             PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,vals2,&B_N,&B_ierr));
6242972d61bSStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
6252972d61bSStefano Zampini             PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,vals2,&B_N,&B_ierr));
6262972d61bSStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
6272972d61bSStefano Zampini             for (j=0;j<B_N;j++) {
6282972d61bSStefano Zampini               for (k=j+1;k<B_N;k++) {
6292972d61bSStefano Zampini                 vals2[k*B_N+j] = vals2[j*B_N+k];
6302972d61bSStefano Zampini               }
6312972d61bSStefano Zampini             }
6322972d61bSStefano Zampini           }
63308122e43SStefano Zampini           ierr = PetscFPTrapPop();CHKERRQ(ierr);
63408122e43SStefano Zampini           ierr = MatSetValues(S_Ej_inv_all,subset_size,dummy_idx,subset_size,dummy_idx,vals2,INSERT_VALUES);CHKERRQ(ierr);
63508122e43SStefano Zampini           ierr = MatDenseRestoreArray(sub_schurs->S_Ej[i],&vals2);CHKERRQ(ierr);
63608122e43SStefano Zampini         }
637a1337663SStefano Zampini         ierr = MatSetValues(S_Ej_tilda_all,subset_size,dummy_idx,subset_size,dummy_idx,vals,INSERT_VALUES);CHKERRQ(ierr);
638b76a6013SStefano Zampini         ierr = MatDenseRestoreArray(Stilda,&vals);CHKERRQ(ierr);
639b76a6013SStefano Zampini         ierr = MatDestroy(&Stilda);CHKERRQ(ierr);
64065d8bf0aSStefano Zampini       }
64165d8bf0aSStefano Zampini       ierr = ISDestroy(&is_E);CHKERRQ(ierr);
642883469d8SStefano Zampini       local_size += subset_size;
643883469d8SStefano Zampini     }
6445ec10c6aSStefano Zampini     ierr = PetscFree2(dummy_idx,work);CHKERRQ(ierr);
6455db18549SStefano Zampini   }
646a1337663SStefano Zampini   ierr = PetscFree(nnz);CHKERRQ(ierr);
647a1337663SStefano Zampini   ierr = MatDestroy(&S_all);CHKERRQ(ierr);
648*d2627357SStefano Zampini   ierr = MatDestroy(&S_all_inv);CHKERRQ(ierr);
6495db18549SStefano Zampini   ierr = MatAssemblyBegin(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
6505db18549SStefano Zampini   ierr = MatAssemblyEnd(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
6518275c6a8SStefano Zampini   if (S_Ej_tilda_all) {
652a1337663SStefano Zampini     ierr = MatAssemblyBegin(S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
653a1337663SStefano Zampini     ierr = MatAssemblyEnd(S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
6548275c6a8SStefano Zampini   }
6558275c6a8SStefano Zampini   if (S_Ej_inv_all) {
65608122e43SStefano Zampini     ierr = MatAssemblyBegin(S_Ej_inv_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
65708122e43SStefano Zampini     ierr = MatAssemblyEnd(S_Ej_inv_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
65808122e43SStefano Zampini   }
659a1337663SStefano Zampini 
6605db18549SStefano Zampini   /* Global matrix of all assembled Schur on subsets */
6615db18549SStefano Zampini   ierr = MatISSetLocalMat(work_mat,sub_schurs->S_Ej_all);CHKERRQ(ierr);
6623927de2eSStefano Zampini   ierr = MatISSetMPIXAIJPreallocation_Private(work_mat,global_schur_subsets,PETSC_TRUE);CHKERRQ(ierr);
6633927de2eSStefano Zampini   ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
6645db18549SStefano Zampini   /* Get local part of (\sum_j S_Ej) */
665883469d8SStefano Zampini   ierr = ISCreateGeneral(PetscObjectComm((PetscObject)sub_schurs->l2gmap),local_size,all_local_idx_G,PETSC_OWN_POINTER,&temp_is);CHKERRQ(ierr);
666d648f858SStefano Zampini   ierr = MatDestroy(&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
667d648f858SStefano Zampini   ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
66808122e43SStefano Zampini 
6698275c6a8SStefano Zampini   if (S_Ej_tilda_all) {
670a1337663SStefano Zampini     ierr = MatISSetLocalMat(work_mat,S_Ej_tilda_all);CHKERRQ(ierr);
671a1337663SStefano Zampini     ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
672d648f858SStefano Zampini     ierr = MatDestroy(&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
673d648f858SStefano Zampini     ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
6748275c6a8SStefano Zampini   }
6758275c6a8SStefano Zampini   if (S_Ej_inv_all) {
67608122e43SStefano Zampini     PetscInt    cum;
67708122e43SStefano Zampini     PetscScalar *array,*array2;
67808122e43SStefano Zampini     ierr = MatISSetLocalMat(work_mat,S_Ej_inv_all);CHKERRQ(ierr);
67908122e43SStefano Zampini     ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
680d648f858SStefano Zampini     ierr = MatDestroy(&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
681d648f858SStefano Zampini     ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
68208122e43SStefano Zampini     /* invert blocks of sum_S_Ej_inv_all */
68308122e43SStefano Zampini     ierr = MatSeqAIJGetArray(sub_schurs->sum_S_Ej_inv_all,&array);CHKERRQ(ierr);
68408122e43SStefano Zampini     ierr = MatSeqAIJGetArray(sub_schurs->sum_S_Ej_tilda_all,&array2);CHKERRQ(ierr);
68508122e43SStefano Zampini     cum = 0;
68608122e43SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
68708122e43SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
68808122e43SStefano Zampini       if (PetscBTLookup(sub_schurs->computed_Stilda_subs,i)) {
68908122e43SStefano Zampini         ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr);
69008122e43SStefano Zampini         ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
6912972d61bSStefano Zampini         if (!sub_schurs->is_hermitian) {
69208122e43SStefano Zampini           PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,array+cum,&B_N,pivots,&B_ierr));
69308122e43SStefano Zampini           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
69406a4b1faSStefano Zampini           PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,array+cum,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
69508122e43SStefano Zampini           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
6962972d61bSStefano Zampini         } else {
6972972d61bSStefano Zampini           PetscInt j,k;
6982972d61bSStefano Zampini 
6992972d61bSStefano Zampini           PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,array+cum,&B_N,&B_ierr));
7002972d61bSStefano Zampini           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
7012972d61bSStefano Zampini           PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,array+cum,&B_N,&B_ierr));
7022972d61bSStefano Zampini           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
7032972d61bSStefano Zampini           for (j=0;j<B_N;j++) {
7042972d61bSStefano Zampini             for (k=j+1;k<B_N;k++) {
7052972d61bSStefano Zampini               array[k*B_N+j+cum] = array[j*B_N+k+cum];
7062972d61bSStefano Zampini             }
7072972d61bSStefano Zampini           }
7082972d61bSStefano Zampini         }
7099552c7c7SStefano Zampini         if (invert_Stildas) { /* Stildas can be singular */
7102972d61bSStefano Zampini           if (!sub_schurs->is_hermitian) {
71108122e43SStefano Zampini             PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,array2+cum,&B_N,pivots,&B_ierr));
71208122e43SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
71306a4b1faSStefano Zampini             PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,array2+cum,&B_N,pivots,Bwork,&B_lwork,&B_ierr));
71408122e43SStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
7152972d61bSStefano Zampini           } else {
7162972d61bSStefano Zampini             PetscInt j,k;
7172972d61bSStefano Zampini 
7182972d61bSStefano Zampini             PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,array2+cum,&B_N,&B_ierr));
7192972d61bSStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
7202972d61bSStefano Zampini             PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,array2+cum,&B_N,&B_ierr));
7212972d61bSStefano Zampini             if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
7222972d61bSStefano Zampini             for (j=0;j<B_N;j++) {
7232972d61bSStefano Zampini               for (k=j+1;k<B_N;k++) {
7242972d61bSStefano Zampini                 array2[k*B_N+j+cum] = array2[j*B_N+k+cum];
7252972d61bSStefano Zampini               }
7262972d61bSStefano Zampini             }
7272972d61bSStefano Zampini           }
7289552c7c7SStefano Zampini         }
72908122e43SStefano Zampini         ierr = PetscFPTrapPop();CHKERRQ(ierr);
73008122e43SStefano Zampini       }
73108122e43SStefano Zampini       cum += subset_size*subset_size;
73208122e43SStefano Zampini     }
73308122e43SStefano Zampini     ierr = MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_inv_all,&array);CHKERRQ(ierr);
73408122e43SStefano Zampini     ierr = MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_tilda_all,&array2);CHKERRQ(ierr);
73508122e43SStefano Zampini   }
7363202ece2SStefano Zampini 
73706a4b1faSStefano Zampini   ierr = PetscFree2(Bwork,pivots);CHKERRQ(ierr);
738a1337663SStefano Zampini   ierr = MatDestroy(&global_schur_subsets);CHKERRQ(ierr);
739a1337663SStefano Zampini   ierr = MatDestroy(&S_Ej_tilda_all);CHKERRQ(ierr);
74008122e43SStefano Zampini   ierr = MatDestroy(&S_Ej_inv_all);CHKERRQ(ierr);
7413202ece2SStefano Zampini   ierr = MatDestroy(&work_mat);CHKERRQ(ierr);
7425db18549SStefano Zampini   ierr = ISDestroy(&temp_is);CHKERRQ(ierr);
743b1b3d7a2SStefano Zampini   PetscFunctionReturn(0);
744b1b3d7a2SStefano Zampini }
745b1b3d7a2SStefano Zampini 
746b1b3d7a2SStefano Zampini #undef __FUNCT__
747b1b3d7a2SStefano Zampini #define __FUNCT__ "PCBDDCSubSchursInit"
7485db18549SStefano Zampini PetscErrorCode PCBDDCSubSchursInit(PCBDDCSubSchurs sub_schurs, Mat A, Mat S, IS is_I, IS is_B, PCBDDCGraph graph, ISLocalToGlobalMapping BtoNmap, PetscInt seqthreshold)
749b1b3d7a2SStefano Zampini {
7509bb4a8caSStefano Zampini   IS                  *faces,*edges,*all_cc,vertices;
751b1b3d7a2SStefano Zampini   PetscInt            *index_sequential,*index_parallel;
752b1b3d7a2SStefano Zampini   PetscInt            *auxlocal_sequential,*auxlocal_parallel;
753b1b3d7a2SStefano Zampini   PetscInt            *auxglobal_sequential,*auxglobal_parallel;
754b96c3477SStefano Zampini   PetscInt            *auxmapping;
755b1b3d7a2SStefano Zampini   PetscInt            i,max_subset_size;
756b1b3d7a2SStefano Zampini   PetscInt            n_sequential_problems,n_local_sequential_problems,n_parallel_problems,n_local_parallel_problems;
757b1b3d7a2SStefano Zampini   PetscInt            n_faces,n_edges,n_all_cc;
758b1b3d7a2SStefano Zampini   PetscBool           is_sorted;
759b1b3d7a2SStefano Zampini   PetscErrorCode  ierr;
760b1b3d7a2SStefano Zampini 
761b1b3d7a2SStefano Zampini   PetscFunctionBegin;
762b1b3d7a2SStefano Zampini   ierr = ISSorted(is_I,&is_sorted);CHKERRQ(ierr);
763b1b3d7a2SStefano Zampini   if (!is_sorted) {
764b1b3d7a2SStefano Zampini     SETERRQ(PetscObjectComm((PetscObject)is_I),PETSC_ERR_PLIB,"IS for I dofs should be shorted");
765b1b3d7a2SStefano Zampini   }
766b1b3d7a2SStefano Zampini   ierr = ISSorted(is_B,&is_sorted);CHKERRQ(ierr);
767b1b3d7a2SStefano Zampini   if (!is_sorted) {
768b1b3d7a2SStefano Zampini     SETERRQ(PetscObjectComm((PetscObject)is_B),PETSC_ERR_PLIB,"IS for B dofs should be shorted");
769b1b3d7a2SStefano Zampini   }
770b1b3d7a2SStefano Zampini 
771b1b3d7a2SStefano Zampini   /* reset any previous data */
772b1b3d7a2SStefano Zampini   ierr = PCBDDCSubSchursReset(sub_schurs);CHKERRQ(ierr);
773b1b3d7a2SStefano Zampini 
774b1b3d7a2SStefano Zampini   /* get index sets for faces and edges */
7759bb4a8caSStefano Zampini   ierr = PCBDDCGraphGetCandidatesIS(graph,&n_faces,&faces,&n_edges,&edges,&vertices);CHKERRQ(ierr);
776b1b3d7a2SStefano Zampini   n_all_cc = n_faces+n_edges;
77708122e43SStefano Zampini   ierr = PetscBTCreate(n_all_cc,&sub_schurs->is_edge);CHKERRQ(ierr);
77808122e43SStefano Zampini   ierr = PetscBTCreate(n_all_cc,&sub_schurs->computed_Stilda_subs);CHKERRQ(ierr);
779b1b3d7a2SStefano Zampini   ierr = PetscMalloc1(n_all_cc,&all_cc);CHKERRQ(ierr);
780b1b3d7a2SStefano Zampini   for (i=0;i<n_faces;i++) {
781b1b3d7a2SStefano Zampini     all_cc[i] = faces[i];
782b1b3d7a2SStefano Zampini   }
783b1b3d7a2SStefano Zampini   for (i=0;i<n_edges;i++) {
784b1b3d7a2SStefano Zampini     all_cc[n_faces+i] = edges[i];
78508122e43SStefano Zampini     ierr = PetscBTSet(sub_schurs->is_edge,n_faces+i);CHKERRQ(ierr);
786b1b3d7a2SStefano Zampini   }
787b1b3d7a2SStefano Zampini   ierr = PetscFree(faces);CHKERRQ(ierr);
788b1b3d7a2SStefano Zampini   ierr = PetscFree(edges);CHKERRQ(ierr);
789b1b3d7a2SStefano Zampini 
790b1b3d7a2SStefano Zampini   /* map interface's subsets */
791b1b3d7a2SStefano Zampini   max_subset_size = 0;
792b1b3d7a2SStefano Zampini   for (i=0;i<n_all_cc;i++) {
793b1b3d7a2SStefano Zampini     PetscInt subset_size;
794b1b3d7a2SStefano Zampini     ierr = ISGetLocalSize(all_cc[i],&subset_size);CHKERRQ(ierr);
795b1b3d7a2SStefano Zampini     max_subset_size = PetscMax(max_subset_size,subset_size);
796b1b3d7a2SStefano Zampini   }
797b1b3d7a2SStefano Zampini   ierr = PetscMalloc1(max_subset_size,&auxmapping);CHKERRQ(ierr);
798b1b3d7a2SStefano Zampini   ierr = PetscMalloc2(graph->ncc,&auxlocal_sequential,
799b1b3d7a2SStefano Zampini                       graph->ncc,&auxlocal_parallel);CHKERRQ(ierr);
800b1b3d7a2SStefano Zampini   ierr = PetscMalloc2(graph->ncc,&index_sequential,
801b1b3d7a2SStefano Zampini                       graph->ncc,&index_parallel);CHKERRQ(ierr);
802b1b3d7a2SStefano Zampini 
803883469d8SStefano Zampini   /* if threshold is negative uses all sequential problems (possibly using MUMPS) */
804883469d8SStefano Zampini   sub_schurs->use_mumps = PETSC_FALSE;
805883469d8SStefano Zampini   if (seqthreshold < 0) {
806883469d8SStefano Zampini     seqthreshold = max_subset_size;
807883469d8SStefano Zampini #if defined(PETSC_HAVE_MUMPS)
8084c6709b3SStefano Zampini     sub_schurs->use_mumps = (PetscBool)(!!A);
809883469d8SStefano Zampini #endif
810883469d8SStefano Zampini   }
811883469d8SStefano Zampini 
812b1b3d7a2SStefano Zampini 
813b1b3d7a2SStefano Zampini   /* determine which problem has to be solved in parallel or sequentially */
814b1b3d7a2SStefano Zampini   n_local_sequential_problems = 0;
815b1b3d7a2SStefano Zampini   n_local_parallel_problems = 0;
816b1b3d7a2SStefano Zampini   for (i=0;i<n_all_cc;i++) {
817b1b3d7a2SStefano Zampini     PetscInt       subset_size,j,min_loc = 0;
818b1b3d7a2SStefano Zampini     const PetscInt *idxs;
819b1b3d7a2SStefano Zampini 
820b1b3d7a2SStefano Zampini     ierr = ISGetLocalSize(all_cc[i],&subset_size);CHKERRQ(ierr);
821b1b3d7a2SStefano Zampini     ierr = ISGetIndices(all_cc[i],&idxs);CHKERRQ(ierr);
822b1b3d7a2SStefano Zampini     ierr = ISLocalToGlobalMappingApply(graph->l2gmap,subset_size,idxs,auxmapping);CHKERRQ(ierr);
823b1b3d7a2SStefano Zampini     for (j=1;j<subset_size;j++) {
824b1b3d7a2SStefano Zampini       if (auxmapping[j]<auxmapping[min_loc]) {
825b1b3d7a2SStefano Zampini         min_loc = j;
826b1b3d7a2SStefano Zampini       }
827b1b3d7a2SStefano Zampini     }
828b1b3d7a2SStefano Zampini     if (subset_size > seqthreshold) {
829b1b3d7a2SStefano Zampini       index_parallel[n_local_parallel_problems] = i;
830b1b3d7a2SStefano Zampini       auxlocal_parallel[n_local_parallel_problems] = idxs[min_loc];
831b1b3d7a2SStefano Zampini       n_local_parallel_problems++;
832b1b3d7a2SStefano Zampini     } else {
833b1b3d7a2SStefano Zampini       index_sequential[n_local_sequential_problems] = i;
834b1b3d7a2SStefano Zampini       auxlocal_sequential[n_local_sequential_problems] = idxs[min_loc];
835b1b3d7a2SStefano Zampini       n_local_sequential_problems++;
836b1b3d7a2SStefano Zampini     }
837b1b3d7a2SStefano Zampini     ierr = ISRestoreIndices(all_cc[i],&idxs);CHKERRQ(ierr);
838b1b3d7a2SStefano Zampini   }
839b1b3d7a2SStefano Zampini 
840b1b3d7a2SStefano Zampini   /* Number parallel problems */
841b1b3d7a2SStefano Zampini   auxglobal_parallel = 0;
842b1b3d7a2SStefano Zampini   ierr = PCBDDCSubsetNumbering(PetscObjectComm((PetscObject)graph->l2gmap),graph->l2gmap,n_local_parallel_problems,auxlocal_parallel,PETSC_NULL,&n_parallel_problems,&auxglobal_parallel);CHKERRQ(ierr);
843b1b3d7a2SStefano Zampini 
844b1b3d7a2SStefano Zampini   /* Number sequential problems */
845b1b3d7a2SStefano Zampini   auxglobal_sequential = 0;
846b1b3d7a2SStefano Zampini   ierr = PCBDDCSubsetNumbering(PetscObjectComm((PetscObject)graph->l2gmap),graph->l2gmap,n_local_sequential_problems,auxlocal_sequential,PETSC_NULL,&n_sequential_problems,&auxglobal_sequential);CHKERRQ(ierr);
847b1b3d7a2SStefano Zampini 
848b1b3d7a2SStefano Zampini   /* update info in sub_schurs */
849aa83b6aeSStefano Zampini   if (A) {
8509ab7bb16SStefano Zampini     PetscBool isseqaij;
8519ab7bb16SStefano Zampini 
8529ab7bb16SStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
8539ab7bb16SStefano Zampini     if (isseqaij) {
8541e9c79c2SStefano Zampini       ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
8551e9c79c2SStefano Zampini       sub_schurs->A = A;
8569ab7bb16SStefano Zampini     } else { /* SeqBAIJ matrices does not support symmetry checking, SeqSBAIJ does not support MatPermute */
8579ab7bb16SStefano Zampini       ierr = MatConvert(A,MATSEQAIJ,MAT_INITIAL_MATRIX,&sub_schurs->A);CHKERRQ(ierr);
8589ab7bb16SStefano Zampini     }
8591e9c79c2SStefano Zampini   }
860b1b3d7a2SStefano Zampini   ierr = PetscObjectReference((PetscObject)S);CHKERRQ(ierr);
861b1b3d7a2SStefano Zampini   sub_schurs->S = S;
862b1b3d7a2SStefano Zampini   ierr = PetscObjectReference((PetscObject)is_I);CHKERRQ(ierr);
863b1b3d7a2SStefano Zampini   sub_schurs->is_I = is_I;
864b1b3d7a2SStefano Zampini   ierr = PetscObjectReference((PetscObject)is_B);CHKERRQ(ierr);
865b1b3d7a2SStefano Zampini   sub_schurs->is_B = is_B;
8665db18549SStefano Zampini   ierr = PetscObjectReference((PetscObject)graph->l2gmap);CHKERRQ(ierr);
8675db18549SStefano Zampini   sub_schurs->l2gmap = graph->l2gmap;
8685db18549SStefano Zampini   ierr = PetscObjectReference((PetscObject)BtoNmap);CHKERRQ(ierr);
8695db18549SStefano Zampini   sub_schurs->BtoNmap = BtoNmap;
870b1b3d7a2SStefano Zampini   sub_schurs->n_subs_seq = n_local_sequential_problems;
871b1b3d7a2SStefano Zampini   sub_schurs->n_subs_par = n_local_parallel_problems;
872b1b3d7a2SStefano Zampini   sub_schurs->n_subs_seq_g = n_sequential_problems;
873b1b3d7a2SStefano Zampini   sub_schurs->n_subs_par_g = n_parallel_problems;
874b1b3d7a2SStefano Zampini   sub_schurs->n_subs = sub_schurs->n_subs_seq + sub_schurs->n_subs_par;
875b1b3d7a2SStefano Zampini   sub_schurs->is_subs = all_cc;
8769bb4a8caSStefano Zampini   if (!sub_schurs->use_mumps) { /* for adaptive selection */
877b96c3477SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
878b96c3477SStefano Zampini       ierr = ISSort(sub_schurs->is_subs[i]);CHKERRQ(ierr);
879b96c3477SStefano Zampini     }
8809bb4a8caSStefano Zampini   }
8819bb4a8caSStefano Zampini   sub_schurs->is_Ej_com = vertices;
882b1b3d7a2SStefano Zampini   sub_schurs->index_sequential = index_sequential;
883b1b3d7a2SStefano Zampini   sub_schurs->index_parallel = index_parallel;
884b1b3d7a2SStefano Zampini   sub_schurs->auxglobal_sequential = auxglobal_sequential;
885b1b3d7a2SStefano Zampini   sub_schurs->auxglobal_parallel = auxglobal_parallel;
886b1b3d7a2SStefano Zampini 
887b96c3477SStefano Zampini 
888b96c3477SStefano Zampini   /* allocate space for schur complements */
889b76a6013SStefano Zampini   ierr = PetscCalloc1(sub_schurs->n_subs,&sub_schurs->S_Ej);CHKERRQ(ierr);
890b96c3477SStefano Zampini   sub_schurs->S_Ej_all = NULL;
891b96c3477SStefano Zampini   sub_schurs->sum_S_Ej_all = NULL;
89208122e43SStefano Zampini   sub_schurs->sum_S_Ej_inv_all = NULL;
893b96c3477SStefano Zampini   sub_schurs->sum_S_Ej_tilda_all = NULL;
894b96c3477SStefano Zampini   sub_schurs->is_Ej_all = NULL;
895b96c3477SStefano Zampini 
896b1b3d7a2SStefano Zampini   /* free workspace */
897b1b3d7a2SStefano Zampini   ierr = PetscFree(auxmapping);CHKERRQ(ierr);
898b1b3d7a2SStefano Zampini   ierr = PetscFree2(auxlocal_sequential,auxlocal_parallel);CHKERRQ(ierr);
899b1b3d7a2SStefano Zampini 
900b1b3d7a2SStefano Zampini   PetscFunctionReturn(0);
901b1b3d7a2SStefano Zampini }
902b1b3d7a2SStefano Zampini 
903b1b3d7a2SStefano Zampini #undef __FUNCT__
90434a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCSubSchursCreate"
90534a97f8cSStefano Zampini PetscErrorCode PCBDDCSubSchursCreate(PCBDDCSubSchurs *sub_schurs)
90634a97f8cSStefano Zampini {
90734a97f8cSStefano Zampini   PCBDDCSubSchurs schurs_ctx;
90834a97f8cSStefano Zampini   PetscErrorCode  ierr;
90934a97f8cSStefano Zampini 
91034a97f8cSStefano Zampini   PetscFunctionBegin;
91134a97f8cSStefano Zampini   ierr = PetscNew(&schurs_ctx);CHKERRQ(ierr);
9125ff63025SStefano Zampini   schurs_ctx->n_subs = 0;
91334a97f8cSStefano Zampini   *sub_schurs = schurs_ctx;
91434a97f8cSStefano Zampini   PetscFunctionReturn(0);
91534a97f8cSStefano Zampini }
91634a97f8cSStefano Zampini 
91734a97f8cSStefano Zampini #undef __FUNCT__
91834a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCSubSchursDestroy"
91934a97f8cSStefano Zampini PetscErrorCode PCBDDCSubSchursDestroy(PCBDDCSubSchurs *sub_schurs)
92034a97f8cSStefano Zampini {
92134a97f8cSStefano Zampini   PetscErrorCode ierr;
92234a97f8cSStefano Zampini 
92334a97f8cSStefano Zampini   PetscFunctionBegin;
92434a97f8cSStefano Zampini   ierr = PCBDDCSubSchursReset(*sub_schurs);CHKERRQ(ierr);
92534a97f8cSStefano Zampini   ierr = PetscFree(*sub_schurs);CHKERRQ(ierr);
92634a97f8cSStefano Zampini   PetscFunctionReturn(0);
92734a97f8cSStefano Zampini }
92834a97f8cSStefano Zampini 
92934a97f8cSStefano Zampini #undef __FUNCT__
93034a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCSubSchursReset"
93134a97f8cSStefano Zampini PetscErrorCode PCBDDCSubSchursReset(PCBDDCSubSchurs sub_schurs)
93234a97f8cSStefano Zampini {
93334a97f8cSStefano Zampini   PetscInt       i;
93434a97f8cSStefano Zampini   PetscErrorCode ierr;
93534a97f8cSStefano Zampini 
93634a97f8cSStefano Zampini   PetscFunctionBegin;
9371e9c79c2SStefano Zampini   ierr = MatDestroy(&sub_schurs->A);CHKERRQ(ierr);
938b1b3d7a2SStefano Zampini   ierr = MatDestroy(&sub_schurs->S);CHKERRQ(ierr);
939b1b3d7a2SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_I);CHKERRQ(ierr);
940b1b3d7a2SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_B);CHKERRQ(ierr);
9415db18549SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->l2gmap);CHKERRQ(ierr);
9425db18549SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->BtoNmap);CHKERRQ(ierr);
94341c3ba1bSStefano Zampini   ierr = MatDestroy(&sub_schurs->S_Ej_all);CHKERRQ(ierr);
94441c3ba1bSStefano Zampini   ierr = MatDestroy(&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
94508122e43SStefano Zampini   ierr = MatDestroy(&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
946a1337663SStefano Zampini   ierr = MatDestroy(&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
9475db18549SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_Ej_all);CHKERRQ(ierr);
9489bb4a8caSStefano Zampini   ierr = ISDestroy(&sub_schurs->is_Ej_com);CHKERRQ(ierr);
94908122e43SStefano Zampini   ierr = PetscBTDestroy(&sub_schurs->is_edge);CHKERRQ(ierr);
95008122e43SStefano Zampini   ierr = PetscBTDestroy(&sub_schurs->computed_Stilda_subs);CHKERRQ(ierr);
95134a97f8cSStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
952b1b3d7a2SStefano Zampini     ierr = ISDestroy(&sub_schurs->is_subs[i]);CHKERRQ(ierr);
95334a97f8cSStefano Zampini     ierr = MatDestroy(&sub_schurs->S_Ej[i]);CHKERRQ(ierr);
95434a97f8cSStefano Zampini   }
95568270318SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_I_layer);CHKERRQ(ierr);
9565ff63025SStefano Zampini   if (sub_schurs->n_subs) {
957b1b3d7a2SStefano Zampini     ierr = PetscFree(sub_schurs->is_subs);CHKERRQ(ierr);
958b76a6013SStefano Zampini     ierr = PetscFree(sub_schurs->S_Ej);CHKERRQ(ierr);
9593dc780c3SStefano Zampini   }
960b1b3d7a2SStefano Zampini   ierr = PetscFree2(sub_schurs->index_sequential,sub_schurs->index_parallel);CHKERRQ(ierr);
961b1b3d7a2SStefano Zampini   ierr = PetscFree(sub_schurs->auxglobal_sequential);CHKERRQ(ierr);
962b1b3d7a2SStefano Zampini   ierr = PetscFree(sub_schurs->auxglobal_parallel);CHKERRQ(ierr);
96334a97f8cSStefano Zampini   sub_schurs->n_subs = 0;
96434a97f8cSStefano Zampini   PetscFunctionReturn(0);
96534a97f8cSStefano Zampini }
96634a97f8cSStefano Zampini 
96734a97f8cSStefano Zampini #undef __FUNCT__
96834a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCAdjGetNextLayer_Private"
9692a155e38SStefano Zampini PETSC_STATIC_INLINE PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt* queue_tip,PetscInt n_prev,PetscBT touched,PetscInt* xadj,PetscInt* adjncy,PetscInt* n_added)
97034a97f8cSStefano Zampini {
97134a97f8cSStefano Zampini   PetscInt       i,j,n;
97234a97f8cSStefano Zampini   PetscErrorCode ierr;
97334a97f8cSStefano Zampini 
97434a97f8cSStefano Zampini   PetscFunctionBegin;
97534a97f8cSStefano Zampini   n = 0;
97634a97f8cSStefano Zampini   for (i=-n_prev;i<0;i++) {
97734a97f8cSStefano Zampini     PetscInt start_dof = queue_tip[i];
97834a97f8cSStefano Zampini     for (j=xadj[start_dof];j<xadj[start_dof+1];j++) {
97934a97f8cSStefano Zampini       PetscInt dof = adjncy[j];
98034a97f8cSStefano Zampini       if (!PetscBTLookup(touched,dof)) {
98134a97f8cSStefano Zampini         ierr = PetscBTSet(touched,dof);CHKERRQ(ierr);
98234a97f8cSStefano Zampini         queue_tip[n] = dof;
98334a97f8cSStefano Zampini         n++;
98434a97f8cSStefano Zampini       }
98534a97f8cSStefano Zampini     }
98634a97f8cSStefano Zampini   }
98734a97f8cSStefano Zampini   *n_added = n;
98834a97f8cSStefano Zampini   PetscFunctionReturn(0);
98934a97f8cSStefano Zampini }
990