xref: /petsc/src/ksp/pc/impls/bddc/bddcschurs.c (revision 883469d828b616a25ec292054b69b8df0fbc2384)
134a97f8cSStefano Zampini #include <../src/ksp/pc/impls/bddc/bddc.h>
234a97f8cSStefano Zampini #include <../src/ksp/pc/impls/bddc/bddcprivate.h>
334a97f8cSStefano Zampini 
434a97f8cSStefano Zampini static PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt*,PetscInt,PetscBT,PetscInt*,PetscInt*,PetscInt*);
534a97f8cSStefano Zampini 
634a97f8cSStefano Zampini #undef __FUNCT__
7b1b3d7a2SStefano Zampini #define __FUNCT__ "PCBDDCSubSchursSetUpNew"
8b1b3d7a2SStefano Zampini PetscErrorCode PCBDDCSubSchursSetUpNew(PCBDDCSubSchurs sub_schurs, PetscInt xadj[], PetscInt adjncy[], PetscInt nlayers)
9b1b3d7a2SStefano Zampini {
10b1b3d7a2SStefano Zampini   Mat                    A_II,A_IB,A_BI,A_BB;
11b1b3d7a2SStefano Zampini   Mat                    AE_II,*AE_IE,*AE_EI,*AE_EE;
12*883469d8SStefano Zampini   Mat                    S_all,global_schur_subsets,work_mat;
135db18549SStefano Zampini   ISLocalToGlobalMapping l2gmap_subsets;
145db18549SStefano Zampini   IS                     is_I,*is_subset_B,temp_is;
15*883469d8SStefano Zampini   PetscInt               *nnz,*all_local_idx_G,*all_local_idx_B,*all_local_idx_N;
165db18549SStefano Zampini   PetscInt               i,subset_size,max_subset_size;
17*883469d8SStefano Zampini   PetscInt               extra,local_size,global_size;
18b1b3d7a2SStefano Zampini   PetscErrorCode         ierr;
19b1b3d7a2SStefano Zampini 
20b1b3d7a2SStefano Zampini   PetscFunctionBegin;
21*883469d8SStefano Zampini 
22b1b3d7a2SStefano Zampini   /* allocate space for schur complements */
2368270318SStefano Zampini   ierr = PetscMalloc4(sub_schurs->n_subs,&sub_schurs->is_AEj_B,
24b1b3d7a2SStefano Zampini                       sub_schurs->n_subs,&sub_schurs->S_Ej,
25b1b3d7a2SStefano Zampini                       sub_schurs->n_subs,&sub_schurs->work1,
26b1b3d7a2SStefano Zampini                       sub_schurs->n_subs,&sub_schurs->work2);CHKERRQ(ierr);
27b1b3d7a2SStefano Zampini 
28b1b3d7a2SStefano Zampini   /* get Schur complement matrices */
29*883469d8SStefano Zampini   if (!sub_schurs->use_mumps) {
30b1b3d7a2SStefano Zampini     ierr = MatSchurComplementGetSubMatrices(sub_schurs->S,&A_II,NULL,&A_IB,&A_BI,&A_BB);CHKERRQ(ierr);
31b1b3d7a2SStefano Zampini     ierr = PetscMalloc4(sub_schurs->n_subs,&is_subset_B,
32b1b3d7a2SStefano Zampini                         sub_schurs->n_subs,&AE_IE,
33b1b3d7a2SStefano Zampini                         sub_schurs->n_subs,&AE_EI,
34b1b3d7a2SStefano Zampini                         sub_schurs->n_subs,&AE_EE);CHKERRQ(ierr);
35b1b3d7a2SStefano Zampini   }
36b1b3d7a2SStefano Zampini 
37b1b3d7a2SStefano Zampini   /* determine interior problems */
38b1b3d7a2SStefano Zampini   if (nlayers >= 0 && xadj != NULL && adjncy != NULL) { /* Interior problems can be different from the original one */
39b1b3d7a2SStefano Zampini     PetscBT                touched;
40b1b3d7a2SStefano Zampini     const PetscInt*        idx_B;
41b1b3d7a2SStefano Zampini     PetscInt               n_I,n_B,n_local_dofs,n_prev_added,j,layer,*local_numbering;
42b1b3d7a2SStefano Zampini 
43b1b3d7a2SStefano Zampini     /* get sizes */
44b1b3d7a2SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr);
45b1b3d7a2SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_B,&n_B);CHKERRQ(ierr);
46b1b3d7a2SStefano Zampini 
47b1b3d7a2SStefano Zampini     ierr = PetscMalloc1(n_I+n_B,&local_numbering);CHKERRQ(ierr);
48b1b3d7a2SStefano Zampini     ierr = PetscBTCreate(n_I+n_B,&touched);CHKERRQ(ierr);
49b1b3d7a2SStefano Zampini     ierr = PetscBTMemzero(n_I+n_B,touched);CHKERRQ(ierr);
50b1b3d7a2SStefano Zampini 
51b1b3d7a2SStefano Zampini     /* all boundary dofs must be skipped when adding layers */
52b1b3d7a2SStefano Zampini     ierr = ISGetIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr);
53b1b3d7a2SStefano Zampini     for (j=0;j<n_B;j++) {
54b1b3d7a2SStefano Zampini       ierr = PetscBTSet(touched,idx_B[j]);CHKERRQ(ierr);
55b1b3d7a2SStefano Zampini     }
56b1b3d7a2SStefano Zampini     ierr = PetscMemcpy(local_numbering,idx_B,n_B*sizeof(PetscInt));CHKERRQ(ierr);
57b1b3d7a2SStefano Zampini     ierr = ISRestoreIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr);
58b1b3d7a2SStefano Zampini 
59b1b3d7a2SStefano Zampini     /* add prescribed number of layers of dofs */
60b1b3d7a2SStefano Zampini     n_local_dofs = n_B;
61b1b3d7a2SStefano Zampini     n_prev_added = n_B;
62b1b3d7a2SStefano Zampini     for (layer=0;layer<nlayers;layer++) {
63b1b3d7a2SStefano Zampini       PetscInt n_added;
64b1b3d7a2SStefano Zampini       if (n_local_dofs == n_I+n_B) break;
65b1b3d7a2SStefano Zampini       if (n_local_dofs > n_I+n_B) {
66b1b3d7a2SStefano 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);
67b1b3d7a2SStefano Zampini       }
68b1b3d7a2SStefano Zampini       ierr = PCBDDCAdjGetNextLayer_Private(local_numbering+n_local_dofs,n_prev_added,touched,xadj,adjncy,&n_added);CHKERRQ(ierr);
69b1b3d7a2SStefano Zampini       n_prev_added = n_added;
70b1b3d7a2SStefano Zampini       n_local_dofs += n_added;
71b1b3d7a2SStefano Zampini       if (!n_added) break;
72b1b3d7a2SStefano Zampini     }
73b1b3d7a2SStefano Zampini     ierr = PetscBTDestroy(&touched);CHKERRQ(ierr);
74b1b3d7a2SStefano Zampini 
75*883469d8SStefano Zampini     /* IS for I layer dofs in original numbering */
7668270318SStefano 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);
77b1b3d7a2SStefano Zampini     ierr = PetscFree(local_numbering);CHKERRQ(ierr);
7868270318SStefano Zampini     ierr = ISSort(sub_schurs->is_I_layer);CHKERRQ(ierr);
79*883469d8SStefano Zampini     /* IS for I layer dofs in I numbering */
80*883469d8SStefano Zampini     if (!sub_schurs->use_mumps) {
81b1b3d7a2SStefano Zampini       ISLocalToGlobalMapping ItoNmap;
82b1b3d7a2SStefano Zampini       ierr = ISLocalToGlobalMappingCreateIS(sub_schurs->is_I,&ItoNmap);CHKERRQ(ierr);
8368270318SStefano Zampini       ierr = ISGlobalToLocalMappingApplyIS(ItoNmap,IS_GTOLM_DROP,sub_schurs->is_I_layer,&is_I);CHKERRQ(ierr);
84b1b3d7a2SStefano Zampini       ierr = ISLocalToGlobalMappingDestroy(&ItoNmap);CHKERRQ(ierr);
85b1b3d7a2SStefano Zampini 
86b1b3d7a2SStefano Zampini       /* II block */
87b1b3d7a2SStefano Zampini       ierr = MatGetSubMatrix(A_II,is_I,is_I,MAT_INITIAL_MATRIX,&AE_II);CHKERRQ(ierr);
88b1b3d7a2SStefano Zampini     }
89b1b3d7a2SStefano Zampini   } else {
90b1b3d7a2SStefano Zampini     PetscInt n_I;
91b1b3d7a2SStefano Zampini 
92b1b3d7a2SStefano Zampini     /* IS for I dofs in original numbering */
93b1b3d7a2SStefano Zampini     ierr = PetscObjectReference((PetscObject)sub_schurs->is_I);CHKERRQ(ierr);
9468270318SStefano Zampini     sub_schurs->is_I_layer = sub_schurs->is_I;
95b1b3d7a2SStefano Zampini 
96b1b3d7a2SStefano Zampini     /* IS for I dofs in I numbering (strided 1) */
97*883469d8SStefano Zampini     if (!sub_schurs->use_mumps) {
98b1b3d7a2SStefano Zampini       ierr = ISGetSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr);
99b1b3d7a2SStefano Zampini       ierr = ISCreateStride(PetscObjectComm((PetscObject)sub_schurs->is_I),n_I,0,1,&is_I);CHKERRQ(ierr);
100b1b3d7a2SStefano Zampini 
101b1b3d7a2SStefano Zampini       /* II block is the same */
102b1b3d7a2SStefano Zampini       ierr = PetscObjectReference((PetscObject)A_II);CHKERRQ(ierr);
103b1b3d7a2SStefano Zampini       AE_II = A_II;
104b1b3d7a2SStefano Zampini     }
105b1b3d7a2SStefano Zampini   }
106b1b3d7a2SStefano Zampini 
107b1b3d7a2SStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
108b1b3d7a2SStefano Zampini     ierr = ISDuplicate(sub_schurs->is_subs[i],&sub_schurs->is_AEj_B[i]);CHKERRQ(ierr);
109b1b3d7a2SStefano Zampini     ierr = ISSort(sub_schurs->is_AEj_B[i]);CHKERRQ(ierr);
110*883469d8SStefano Zampini   }
111*883469d8SStefano Zampini 
112*883469d8SStefano Zampini   /* Get info on subset sizes and sum of all subsets sizes */
113*883469d8SStefano Zampini   max_subset_size = 0;
114*883469d8SStefano Zampini   local_size = 0;
115*883469d8SStefano Zampini   for (i=0;i<sub_schurs->n_subs_seq;i++) {
116*883469d8SStefano Zampini     PetscInt j = sub_schurs->index_sequential[i];
117*883469d8SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_AEj_B[j],&subset_size);CHKERRQ(ierr);
118*883469d8SStefano Zampini     max_subset_size = PetscMax(subset_size,max_subset_size);
119*883469d8SStefano Zampini     local_size += subset_size;
120*883469d8SStefano Zampini   }
121*883469d8SStefano Zampini 
122*883469d8SStefano Zampini   /* Work arrays for local indices */
123*883469d8SStefano Zampini   ierr = PetscMalloc1(local_size,&all_local_idx_B);CHKERRQ(ierr);
124*883469d8SStefano Zampini   extra = 0;
125*883469d8SStefano Zampini   if (sub_schurs->use_mumps) {
126*883469d8SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_I_layer,&extra);CHKERRQ(ierr);
127*883469d8SStefano Zampini   }
128*883469d8SStefano Zampini   ierr = PetscMalloc1(local_size+extra,&all_local_idx_N);CHKERRQ(ierr);
129*883469d8SStefano Zampini   if (extra) {
130*883469d8SStefano Zampini     const PetscInt *idxs;
131*883469d8SStefano Zampini     ierr = ISGetIndices(sub_schurs->is_I_layer,&idxs);CHKERRQ(ierr);
132*883469d8SStefano Zampini     ierr = PetscMemcpy(all_local_idx_N,idxs,extra*sizeof(PetscInt));CHKERRQ(ierr);
133*883469d8SStefano Zampini     ierr = ISRestoreIndices(sub_schurs->is_I_layer,&idxs);CHKERRQ(ierr);
134*883469d8SStefano Zampini   }
135*883469d8SStefano Zampini   ierr = PetscMalloc1(local_size,&nnz);CHKERRQ(ierr);
136*883469d8SStefano Zampini 
137*883469d8SStefano Zampini   /* Get local indices in local numbering */
138*883469d8SStefano Zampini   local_size = 0;
139*883469d8SStefano Zampini   for (i=0;i<sub_schurs->n_subs_seq;i++) {
140*883469d8SStefano Zampini     PetscInt j;
141*883469d8SStefano Zampini     const    PetscInt *idxs;
142*883469d8SStefano Zampini 
143*883469d8SStefano Zampini     PetscInt local_problem_index = sub_schurs->index_sequential[i];
144*883469d8SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_AEj_B[local_problem_index],&subset_size);CHKERRQ(ierr);
145*883469d8SStefano Zampini     ierr = ISGetIndices(sub_schurs->is_AEj_B[local_problem_index],&idxs);CHKERRQ(ierr);
146*883469d8SStefano Zampini     /* subset indices in local numbering */
147*883469d8SStefano Zampini     ierr = PetscMemcpy(all_local_idx_N+local_size+extra,idxs,subset_size*sizeof(PetscInt));CHKERRQ(ierr);
148*883469d8SStefano Zampini     ierr = ISRestoreIndices(sub_schurs->is_AEj_B[local_problem_index],&idxs);CHKERRQ(ierr);
149*883469d8SStefano Zampini     for (j=0;j<subset_size;j++) nnz[local_size+j] = subset_size;
150*883469d8SStefano Zampini     local_size += subset_size;
151*883469d8SStefano Zampini   }
152*883469d8SStefano Zampini 
153*883469d8SStefano Zampini   S_all = NULL;
154*883469d8SStefano Zampini   if (!sub_schurs->use_mumps) {
155*883469d8SStefano Zampini     /* subsets in original and boundary numbering */
156*883469d8SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
1575db18549SStefano Zampini       ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,sub_schurs->is_AEj_B[i],&is_subset_B[i]);CHKERRQ(ierr);
158b1b3d7a2SStefano Zampini     }
159b1b3d7a2SStefano Zampini 
160b1b3d7a2SStefano Zampini     /* EE block */
161b1b3d7a2SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
162b1b3d7a2SStefano Zampini       ierr = MatGetSubMatrix(A_BB,is_subset_B[i],is_subset_B[i],MAT_INITIAL_MATRIX,&AE_EE[i]);CHKERRQ(ierr);
163b1b3d7a2SStefano Zampini     }
164b1b3d7a2SStefano Zampini     /* IE block */
165b1b3d7a2SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
166b1b3d7a2SStefano Zampini       ierr = MatGetSubMatrix(A_IB,is_I,is_subset_B[i],MAT_INITIAL_MATRIX,&AE_IE[i]);CHKERRQ(ierr);
167b1b3d7a2SStefano Zampini     }
168b1b3d7a2SStefano Zampini     /* EI block */
169b1b3d7a2SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
170b1b3d7a2SStefano Zampini       ierr = MatGetSubMatrix(A_BI,is_subset_B[i],is_I,MAT_INITIAL_MATRIX,&AE_EI[i]);CHKERRQ(ierr);
171b1b3d7a2SStefano Zampini     }
172b1b3d7a2SStefano Zampini 
173b1b3d7a2SStefano Zampini     /* setup Schur complements on subset */
174b1b3d7a2SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
175b1b3d7a2SStefano Zampini       ierr = MatCreateSchurComplement(AE_II,AE_II,AE_IE[i],AE_EI[i],AE_EE[i],&sub_schurs->S_Ej[i]);CHKERRQ(ierr);
176b1b3d7a2SStefano Zampini       if (AE_II == A_II) { /* we can reuse the same ksp */
177b1b3d7a2SStefano Zampini         KSP ksp;
178b1b3d7a2SStefano Zampini         ierr = MatSchurComplementGetKSP(sub_schurs->S,&ksp);CHKERRQ(ierr);
179b1b3d7a2SStefano Zampini         ierr = MatSchurComplementSetKSP(sub_schurs->S_Ej[i],ksp);CHKERRQ(ierr);
180b1b3d7a2SStefano Zampini       } else { /* build new ksp object which inherits ksp and pc types from the original one */
181b1b3d7a2SStefano Zampini         KSP      origksp,schurksp;
182b1b3d7a2SStefano Zampini         PC       origpc,schurpc;
183b1b3d7a2SStefano Zampini         KSPType  ksp_type;
184b1b3d7a2SStefano Zampini         PCType   pc_type;
185b1b3d7a2SStefano Zampini         PetscInt n_internal;
186b1b3d7a2SStefano Zampini 
187b1b3d7a2SStefano Zampini         ierr = MatSchurComplementGetKSP(sub_schurs->S,&origksp);CHKERRQ(ierr);
188b1b3d7a2SStefano Zampini         ierr = MatSchurComplementGetKSP(sub_schurs->S_Ej[i],&schurksp);CHKERRQ(ierr);
189b1b3d7a2SStefano Zampini         ierr = KSPGetType(origksp,&ksp_type);CHKERRQ(ierr);
190b1b3d7a2SStefano Zampini         ierr = KSPSetType(schurksp,ksp_type);CHKERRQ(ierr);
191b1b3d7a2SStefano Zampini         ierr = KSPGetPC(schurksp,&schurpc);CHKERRQ(ierr);
192b1b3d7a2SStefano Zampini         ierr = KSPGetPC(origksp,&origpc);CHKERRQ(ierr);
193b1b3d7a2SStefano Zampini         ierr = PCGetType(origpc,&pc_type);CHKERRQ(ierr);
194b1b3d7a2SStefano Zampini         ierr = PCSetType(schurpc,pc_type);CHKERRQ(ierr);
195b1b3d7a2SStefano Zampini         ierr = ISGetSize(is_I,&n_internal);CHKERRQ(ierr);
196b1b3d7a2SStefano Zampini         if (n_internal) { /* UMFPACK gives error with 0 sized problems */
197b1b3d7a2SStefano Zampini           MatSolverPackage solver=NULL;
198b1b3d7a2SStefano Zampini           ierr = PCFactorGetMatSolverPackage(origpc,(const MatSolverPackage*)&solver);CHKERRQ(ierr);
199b1b3d7a2SStefano Zampini           if (solver) {
200b1b3d7a2SStefano Zampini             ierr = PCFactorSetMatSolverPackage(schurpc,solver);CHKERRQ(ierr);
201b1b3d7a2SStefano Zampini           }
202b1b3d7a2SStefano Zampini         }
203b1b3d7a2SStefano Zampini         ierr = KSPSetUp(schurksp);CHKERRQ(ierr);
204b1b3d7a2SStefano Zampini       }
205b1b3d7a2SStefano Zampini     }
206b1b3d7a2SStefano Zampini     /* free */
207b1b3d7a2SStefano Zampini     ierr = ISDestroy(&is_I);CHKERRQ(ierr);
208b1b3d7a2SStefano Zampini     ierr = MatDestroy(&AE_II);CHKERRQ(ierr);
209b1b3d7a2SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
210b1b3d7a2SStefano Zampini       ierr = MatDestroy(&AE_EE[i]);CHKERRQ(ierr);
211b1b3d7a2SStefano Zampini       ierr = MatDestroy(&AE_IE[i]);CHKERRQ(ierr);
212b1b3d7a2SStefano Zampini       ierr = MatDestroy(&AE_EI[i]);CHKERRQ(ierr);
213b1b3d7a2SStefano Zampini       ierr = ISDestroy(&is_subset_B[i]);CHKERRQ(ierr);
214b1b3d7a2SStefano Zampini     }
215b1b3d7a2SStefano Zampini     ierr = PetscFree4(is_subset_B,AE_IE,AE_EI,AE_EE);CHKERRQ(ierr);
216*883469d8SStefano Zampini #if defined(PETSC_HAVE_MUMPS)
217*883469d8SStefano Zampini   } else {
218*883469d8SStefano Zampini     Mat           A,F;
219*883469d8SStefano Zampini     IS            is_A_all;
220*883469d8SStefano Zampini     PetscBool     is_symmetric;
221*883469d8SStefano Zampini     PetscInt      *idxs_schur,n_I;
222*883469d8SStefano Zampini 
223*883469d8SStefano Zampini     /* get working mat */
224*883469d8SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_I_layer,&n_I);CHKERRQ(ierr);
225*883469d8SStefano Zampini     ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size+n_I,all_local_idx_N,PETSC_COPY_VALUES,&is_A_all);CHKERRQ(ierr);
226*883469d8SStefano Zampini     ierr = MatGetSubMatrixUnsorted(sub_schurs->A,is_A_all,is_A_all,&A);CHKERRQ(ierr);
227*883469d8SStefano Zampini     ierr = ISDestroy(&is_A_all);CHKERRQ(ierr);
228*883469d8SStefano Zampini 
229*883469d8SStefano Zampini     ierr = MatIsSymmetric(sub_schurs->A,0.0,&is_symmetric);CHKERRQ(ierr);
230*883469d8SStefano Zampini     if (is_symmetric) {
231*883469d8SStefano Zampini       ierr = MatGetFactor(A,MATSOLVERMUMPS,MAT_FACTOR_CHOLESKY,&F);CHKERRQ(ierr);
232*883469d8SStefano Zampini     } else {
233*883469d8SStefano Zampini       ierr = MatGetFactor(A,MATSOLVERMUMPS,MAT_FACTOR_LU,&F);CHKERRQ(ierr);
234*883469d8SStefano Zampini     }
235*883469d8SStefano Zampini 
236*883469d8SStefano Zampini     /* subsets ordered last */
237*883469d8SStefano Zampini     ierr = PetscMalloc1(local_size,&idxs_schur);CHKERRQ(ierr);
238*883469d8SStefano Zampini     for (i=0;i<local_size;i++) {
239*883469d8SStefano Zampini       idxs_schur[i] = n_I+i+1;
240*883469d8SStefano Zampini     }
241*883469d8SStefano Zampini     ierr = MatMumpsSetSchurIndices(F,local_size,idxs_schur);CHKERRQ(ierr);
242*883469d8SStefano Zampini     ierr = PetscFree(idxs_schur);CHKERRQ(ierr);
243*883469d8SStefano Zampini 
244*883469d8SStefano Zampini     /* factorization step */
245*883469d8SStefano Zampini     if (is_symmetric) {
246*883469d8SStefano Zampini       ierr = MatCholeskyFactorSymbolic(F,A,NULL,NULL);CHKERRQ(ierr);
247*883469d8SStefano Zampini       ierr = MatCholeskyFactorNumeric(F,A,NULL);CHKERRQ(ierr);
248*883469d8SStefano Zampini     } else {
249*883469d8SStefano Zampini       ierr = MatLUFactorSymbolic(F,A,NULL,NULL,NULL);CHKERRQ(ierr);
250*883469d8SStefano Zampini       ierr = MatLUFactorNumeric(F,A,NULL);CHKERRQ(ierr);
251*883469d8SStefano Zampini     }
252*883469d8SStefano Zampini 
253*883469d8SStefano Zampini     /* get explicit Schur Complement computed during numeric factorization */
254*883469d8SStefano Zampini     ierr = MatMumpsGetSchurComplement(F,&S_all);CHKERRQ(ierr);
255*883469d8SStefano Zampini 
256*883469d8SStefano Zampini     /* free workspace */
257*883469d8SStefano Zampini     ierr = MatDestroy(&A);CHKERRQ(ierr);
258*883469d8SStefano Zampini     ierr = MatDestroy(&F);CHKERRQ(ierr);
259*883469d8SStefano Zampini 
260*883469d8SStefano Zampini     /* unused */
261*883469d8SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
262*883469d8SStefano Zampini        sub_schurs->S_Ej[i] = 0;
263*883469d8SStefano Zampini     }
264*883469d8SStefano Zampini #endif
265b1b3d7a2SStefano Zampini   }
2665db18549SStefano Zampini 
2679221af80SStefano Zampini   /* TODO: just for compatibility with the previous version, needs to be fixed */
2689221af80SStefano Zampini   for (i=0;i<sub_schurs->n_subs_par;i++) {
2699221af80SStefano Zampini     PetscInt j = sub_schurs->index_parallel[i];
2709221af80SStefano Zampini     ierr = MatCreateVecs(sub_schurs->S_Ej[j],&sub_schurs->work1[j],&sub_schurs->work2[j]);CHKERRQ(ierr);
2719221af80SStefano Zampini   }
2729221af80SStefano Zampini   for (i=0;i<sub_schurs->n_subs_seq;i++) {
2739221af80SStefano Zampini      sub_schurs->work1[sub_schurs->index_sequential[i]] = 0;
2749221af80SStefano Zampini      sub_schurs->work2[sub_schurs->index_sequential[i]] = 0;
2759221af80SStefano Zampini   }
2769221af80SStefano Zampini 
2775db18549SStefano Zampini   if (!sub_schurs->n_subs_seq_g) {
2789221af80SStefano Zampini     sub_schurs->S_Ej_all = 0;
2799221af80SStefano Zampini     sub_schurs->sum_S_Ej_all = 0;
2809221af80SStefano Zampini     sub_schurs->is_Ej_all = 0;
2815db18549SStefano Zampini     PetscFunctionReturn(0);
2825db18549SStefano Zampini   }
2835db18549SStefano Zampini 
2845db18549SStefano Zampini   /* subset indices in local boundary numbering */
285*883469d8SStefano Zampini   ierr = ISGlobalToLocalMappingApply(sub_schurs->BtoNmap,IS_GTOLM_DROP,local_size,all_local_idx_N+extra,&subset_size,all_local_idx_B);CHKERRQ(ierr);
2865db18549SStefano Zampini   if (subset_size != local_size) {
2875db18549SStefano Zampini     SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in sub_schurs serial (BtoNmap)! %d != %d\n",subset_size,local_size);
2885db18549SStefano Zampini   }
2895db18549SStefano Zampini   ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size,all_local_idx_B,PETSC_OWN_POINTER,&sub_schurs->is_Ej_all);CHKERRQ(ierr);
2905db18549SStefano Zampini 
2915db18549SStefano Zampini   /* Local matrix of all local Schur on subsets */
2925db18549SStefano Zampini   ierr = MatCreate(PETSC_COMM_SELF,&sub_schurs->S_Ej_all);CHKERRQ(ierr);
2935db18549SStefano Zampini   ierr = MatSetSizes(sub_schurs->S_Ej_all,PETSC_DECIDE,PETSC_DECIDE,local_size,local_size);CHKERRQ(ierr);
2945db18549SStefano Zampini   ierr = MatSetType(sub_schurs->S_Ej_all,MATAIJ);CHKERRQ(ierr);
2955db18549SStefano Zampini   ierr = MatSeqAIJSetPreallocation(sub_schurs->S_Ej_all,0,nnz);CHKERRQ(ierr);
2965db18549SStefano Zampini   ierr = PetscFree(nnz);CHKERRQ(ierr);
2975db18549SStefano Zampini 
298*883469d8SStefano Zampini   if (!sub_schurs->use_mumps) {
2995db18549SStefano Zampini     PetscScalar *fill_vals;
3005db18549SStefano Zampini     PetscInt    *dummy_idx;
3015db18549SStefano Zampini 
3025db18549SStefano Zampini     /* Work arrays */
3035db18549SStefano Zampini     ierr = PetscMalloc2(max_subset_size,&dummy_idx,max_subset_size*max_subset_size,&fill_vals);CHKERRQ(ierr);
3045db18549SStefano Zampini 
3055db18549SStefano Zampini     /* Loop on local problems to compute Schur complements explicitly (TODO; optimize)*/
3065db18549SStefano Zampini     local_size = 0;
3075db18549SStefano Zampini     for (i=0;i<sub_schurs->n_subs_seq;i++) {
3085db18549SStefano Zampini       Vec work1,work2;
3095db18549SStefano Zampini       PetscInt j,local_problem_index = sub_schurs->index_sequential[i];
3105db18549SStefano Zampini 
3119221af80SStefano Zampini       ierr = MatCreateVecs(sub_schurs->S_Ej[local_problem_index],&work1,&work2);CHKERRQ(ierr);
3125db18549SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_AEj_B[local_problem_index],&subset_size);CHKERRQ(ierr);
3135db18549SStefano Zampini       /* local Schur */
3145db18549SStefano Zampini       for (j=0;j<subset_size;j++) {
3155db18549SStefano Zampini         ierr = VecSet(work1,0.0);CHKERRQ(ierr);
3165db18549SStefano Zampini         ierr = VecSetValue(work1,j,1.0,INSERT_VALUES);CHKERRQ(ierr);
3175db18549SStefano Zampini         ierr = VecPlaceArray(work2,&fill_vals[j*subset_size]);CHKERRQ(ierr);
3185db18549SStefano Zampini         ierr = MatMult(sub_schurs->S_Ej[local_problem_index],work1,work2);CHKERRQ(ierr);
3195db18549SStefano Zampini         ierr = VecResetArray(work2);CHKERRQ(ierr);
3205db18549SStefano Zampini       }
3215db18549SStefano Zampini       for (j=0;j<subset_size;j++) {
3225db18549SStefano Zampini         dummy_idx[j]=local_size+j;
3235db18549SStefano Zampini       }
3245db18549SStefano Zampini       ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,fill_vals,INSERT_VALUES);CHKERRQ(ierr);
3255db18549SStefano Zampini       local_size += subset_size;
3265db18549SStefano Zampini       ierr = VecDestroy(&work1);CHKERRQ(ierr);
3275db18549SStefano Zampini       ierr = VecDestroy(&work2);CHKERRQ(ierr);
3285db18549SStefano Zampini     }
3295db18549SStefano Zampini     ierr = PetscFree2(dummy_idx,fill_vals);CHKERRQ(ierr);
330*883469d8SStefano Zampini   } else {
331*883469d8SStefano Zampini     PetscInt    *dummy_idx,n_all;
332*883469d8SStefano Zampini     PetscScalar *vals,*fill_vals;
333*883469d8SStefano Zampini 
334*883469d8SStefano Zampini     /* Work arrays */
335*883469d8SStefano Zampini     ierr = PetscMalloc1(max_subset_size,&dummy_idx);CHKERRQ(ierr);
336*883469d8SStefano Zampini     ierr = MatGetSize(S_all,&n_all,NULL);CHKERRQ(ierr);
337*883469d8SStefano Zampini     ierr = MatDenseGetArray(S_all,&vals);CHKERRQ(ierr);
338*883469d8SStefano Zampini     local_size = 0;
339*883469d8SStefano Zampini     subset_size = 0;
340*883469d8SStefano Zampini     fill_vals = vals;
341*883469d8SStefano Zampini     for (i=0;i<sub_schurs->n_subs_seq;i++) {
342*883469d8SStefano Zampini       PetscInt j,lpi = sub_schurs->index_sequential[i];
343*883469d8SStefano Zampini 
344*883469d8SStefano Zampini       fill_vals += subset_size;
345*883469d8SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_AEj_B[lpi],&subset_size);CHKERRQ(ierr);
346*883469d8SStefano Zampini       for (j=0;j<subset_size;j++) {
347*883469d8SStefano Zampini         dummy_idx[j]=local_size+j;
348*883469d8SStefano Zampini       }
349*883469d8SStefano Zampini       for (j=0;j<subset_size;j++) {
350*883469d8SStefano Zampini         ierr = MatSetValues(sub_schurs->S_Ej_all,1,dummy_idx+j,subset_size,dummy_idx,fill_vals,INSERT_VALUES);CHKERRQ(ierr);
351*883469d8SStefano Zampini         fill_vals += n_all;
352*883469d8SStefano Zampini       }
353*883469d8SStefano Zampini       local_size += subset_size;
354*883469d8SStefano Zampini     }
355*883469d8SStefano Zampini     ierr = MatDenseRestoreArray(S_all,&vals);CHKERRQ(ierr);
356*883469d8SStefano Zampini     ierr = PetscFree(dummy_idx);CHKERRQ(ierr);
3575db18549SStefano Zampini   }
3585db18549SStefano Zampini   ierr = MatAssemblyBegin(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
3595db18549SStefano Zampini   ierr = MatAssemblyEnd(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
360*883469d8SStefano Zampini   ierr = MatDestroy(&S_all);CHKERRQ(ierr);
3615db18549SStefano Zampini 
3625db18549SStefano Zampini   /* Global matrix of all assembled Schur on subsets */
363*883469d8SStefano 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);
3645db18549SStefano Zampini   ierr = ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)sub_schurs->l2gmap),1,local_size,all_local_idx_G,PETSC_COPY_VALUES,&l2gmap_subsets);CHKERRQ(ierr);
3655db18549SStefano Zampini   ierr = MatCreateIS(PetscObjectComm((PetscObject)sub_schurs->l2gmap),1,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size,l2gmap_subsets,&work_mat);CHKERRQ(ierr);
3665db18549SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&l2gmap_subsets);CHKERRQ(ierr);
3675db18549SStefano Zampini   ierr = MatISSetLocalMat(work_mat,sub_schurs->S_Ej_all);CHKERRQ(ierr);
3685db18549SStefano Zampini   ierr = MatISGetMPIXAIJ(work_mat,MAT_INITIAL_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
3695db18549SStefano Zampini   ierr = MatDestroy(&work_mat);CHKERRQ(ierr);
3705db18549SStefano Zampini 
3715db18549SStefano Zampini   /* Get local part of (\sum_j S_Ej) */
372*883469d8SStefano Zampini   ierr = PetscFree(all_local_idx_N);CHKERRQ(ierr);
373*883469d8SStefano Zampini   ierr = ISCreateGeneral(PetscObjectComm((PetscObject)sub_schurs->l2gmap),local_size,all_local_idx_G,PETSC_OWN_POINTER,&temp_is);CHKERRQ(ierr);
374*883469d8SStefano Zampini   ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
3755db18549SStefano Zampini   ierr = MatDestroy(&global_schur_subsets);CHKERRQ(ierr);
3765db18549SStefano Zampini   ierr = ISDestroy(&temp_is);CHKERRQ(ierr);
377b1b3d7a2SStefano Zampini   PetscFunctionReturn(0);
378b1b3d7a2SStefano Zampini }
379b1b3d7a2SStefano Zampini 
380b1b3d7a2SStefano Zampini #undef __FUNCT__
381b1b3d7a2SStefano Zampini #define __FUNCT__ "PCBDDCSubSchursInit"
3825db18549SStefano Zampini PetscErrorCode PCBDDCSubSchursInit(PCBDDCSubSchurs sub_schurs, Mat A, Mat S, IS is_I, IS is_B, PCBDDCGraph graph, ISLocalToGlobalMapping BtoNmap, PetscInt seqthreshold)
383b1b3d7a2SStefano Zampini {
384b1b3d7a2SStefano Zampini   IS                  *faces,*edges,*all_cc;
385b1b3d7a2SStefano Zampini   PetscInt            *index_sequential,*index_parallel;
386b1b3d7a2SStefano Zampini   PetscInt            *auxlocal_sequential,*auxlocal_parallel;
387b1b3d7a2SStefano Zampini   PetscInt            *auxglobal_sequential,*auxglobal_parallel;
388b1b3d7a2SStefano Zampini   PetscInt            *auxmapping;//,*idxs;
389b1b3d7a2SStefano Zampini   PetscInt            i,max_subset_size;
390b1b3d7a2SStefano Zampini   PetscInt            n_sequential_problems,n_local_sequential_problems,n_parallel_problems,n_local_parallel_problems;
391b1b3d7a2SStefano Zampini   PetscInt            n_faces,n_edges,n_all_cc;
392b1b3d7a2SStefano Zampini   PetscBool           is_sorted;
393b1b3d7a2SStefano Zampini   PetscErrorCode  ierr;
394b1b3d7a2SStefano Zampini 
395b1b3d7a2SStefano Zampini   PetscFunctionBegin;
396b1b3d7a2SStefano Zampini   ierr = ISSorted(is_I,&is_sorted);CHKERRQ(ierr);
397b1b3d7a2SStefano Zampini   if (!is_sorted) {
398b1b3d7a2SStefano Zampini     SETERRQ(PetscObjectComm((PetscObject)is_I),PETSC_ERR_PLIB,"IS for I dofs should be shorted");
399b1b3d7a2SStefano Zampini   }
400b1b3d7a2SStefano Zampini   ierr = ISSorted(is_B,&is_sorted);CHKERRQ(ierr);
401b1b3d7a2SStefano Zampini   if (!is_sorted) {
402b1b3d7a2SStefano Zampini     SETERRQ(PetscObjectComm((PetscObject)is_B),PETSC_ERR_PLIB,"IS for B dofs should be shorted");
403b1b3d7a2SStefano Zampini   }
404b1b3d7a2SStefano Zampini 
405b1b3d7a2SStefano Zampini   /* reset any previous data */
406b1b3d7a2SStefano Zampini   ierr = PCBDDCSubSchursReset(sub_schurs);CHKERRQ(ierr);
407b1b3d7a2SStefano Zampini 
408b1b3d7a2SStefano Zampini   /* get index sets for faces and edges */
409b1b3d7a2SStefano Zampini   ierr = PCBDDCGraphGetCandidatesIS(graph,&n_faces,&faces,&n_edges,&edges,NULL);CHKERRQ(ierr);
410b1b3d7a2SStefano Zampini   n_all_cc = n_faces+n_edges;
411b1b3d7a2SStefano Zampini   ierr = PetscMalloc1(n_all_cc,&all_cc);CHKERRQ(ierr);
412b1b3d7a2SStefano Zampini   for (i=0;i<n_faces;i++) {
413b1b3d7a2SStefano Zampini     all_cc[i] = faces[i];
414b1b3d7a2SStefano Zampini   }
415b1b3d7a2SStefano Zampini   for (i=0;i<n_edges;i++) {
416b1b3d7a2SStefano Zampini     all_cc[n_faces+i] = edges[i];
417b1b3d7a2SStefano Zampini   }
418b1b3d7a2SStefano Zampini   ierr = PetscFree(faces);CHKERRQ(ierr);
419b1b3d7a2SStefano Zampini   ierr = PetscFree(edges);CHKERRQ(ierr);
420b1b3d7a2SStefano Zampini 
421b1b3d7a2SStefano Zampini   /* map interface's subsets */
422b1b3d7a2SStefano Zampini   max_subset_size = 0;
423b1b3d7a2SStefano Zampini   for (i=0;i<n_all_cc;i++) {
424b1b3d7a2SStefano Zampini     PetscInt subset_size;
425b1b3d7a2SStefano Zampini     ierr = ISGetLocalSize(all_cc[i],&subset_size);CHKERRQ(ierr);
426b1b3d7a2SStefano Zampini     max_subset_size = PetscMax(max_subset_size,subset_size);
427b1b3d7a2SStefano Zampini   }
428b1b3d7a2SStefano Zampini   ierr = PetscMalloc1(max_subset_size,&auxmapping);CHKERRQ(ierr);
429b1b3d7a2SStefano Zampini   ierr = PetscMalloc2(graph->ncc,&auxlocal_sequential,
430b1b3d7a2SStefano Zampini                       graph->ncc,&auxlocal_parallel);CHKERRQ(ierr);
431b1b3d7a2SStefano Zampini   ierr = PetscMalloc2(graph->ncc,&index_sequential,
432b1b3d7a2SStefano Zampini                       graph->ncc,&index_parallel);CHKERRQ(ierr);
433b1b3d7a2SStefano Zampini 
434*883469d8SStefano Zampini   /* if threshold is negative uses all sequential problems (possibly using MUMPS) */
435*883469d8SStefano Zampini   sub_schurs->use_mumps = PETSC_FALSE;
436*883469d8SStefano Zampini   if (seqthreshold < 0) {
437*883469d8SStefano Zampini     seqthreshold = max_subset_size;
438*883469d8SStefano Zampini #if defined(PETSC_HAVE_MUMPS)
439*883469d8SStefano Zampini     sub_schurs->use_mumps = !!A;
440*883469d8SStefano Zampini #endif
441*883469d8SStefano Zampini   }
442*883469d8SStefano Zampini 
443b1b3d7a2SStefano Zampini 
444b1b3d7a2SStefano Zampini   /* determine which problem has to be solved in parallel or sequentially */
445b1b3d7a2SStefano Zampini   n_local_sequential_problems = 0;
446b1b3d7a2SStefano Zampini   n_local_parallel_problems = 0;
447b1b3d7a2SStefano Zampini   for (i=0;i<n_all_cc;i++) {
448b1b3d7a2SStefano Zampini     PetscInt       subset_size,j,min_loc = 0;
449b1b3d7a2SStefano Zampini     const PetscInt *idxs;
450b1b3d7a2SStefano Zampini 
451b1b3d7a2SStefano Zampini     ierr = ISGetLocalSize(all_cc[i],&subset_size);CHKERRQ(ierr);
452b1b3d7a2SStefano Zampini     ierr = ISGetIndices(all_cc[i],&idxs);CHKERRQ(ierr);
453b1b3d7a2SStefano Zampini     ierr = ISLocalToGlobalMappingApply(graph->l2gmap,subset_size,idxs,auxmapping);CHKERRQ(ierr);
454b1b3d7a2SStefano Zampini     for (j=1;j<subset_size;j++) {
455b1b3d7a2SStefano Zampini       if (auxmapping[j]<auxmapping[min_loc]) {
456b1b3d7a2SStefano Zampini         min_loc = j;
457b1b3d7a2SStefano Zampini       }
458b1b3d7a2SStefano Zampini     }
459b1b3d7a2SStefano Zampini     if (subset_size > seqthreshold) {
460b1b3d7a2SStefano Zampini       index_parallel[n_local_parallel_problems] = i;
461b1b3d7a2SStefano Zampini       auxlocal_parallel[n_local_parallel_problems] = idxs[min_loc];
462b1b3d7a2SStefano Zampini       n_local_parallel_problems++;
463b1b3d7a2SStefano Zampini     } else {
464b1b3d7a2SStefano Zampini       index_sequential[n_local_sequential_problems] = i;
465b1b3d7a2SStefano Zampini       auxlocal_sequential[n_local_sequential_problems] = idxs[min_loc];
466b1b3d7a2SStefano Zampini       n_local_sequential_problems++;
467b1b3d7a2SStefano Zampini     }
468b1b3d7a2SStefano Zampini     ierr = ISRestoreIndices(all_cc[i],&idxs);CHKERRQ(ierr);
469b1b3d7a2SStefano Zampini   }
470b1b3d7a2SStefano Zampini 
471b1b3d7a2SStefano Zampini   /* Number parallel problems */
472b1b3d7a2SStefano Zampini   auxglobal_parallel = 0;
473b1b3d7a2SStefano Zampini   ierr = PCBDDCSubsetNumbering(PetscObjectComm((PetscObject)graph->l2gmap),graph->l2gmap,n_local_parallel_problems,auxlocal_parallel,PETSC_NULL,&n_parallel_problems,&auxglobal_parallel);CHKERRQ(ierr);
474b1b3d7a2SStefano Zampini 
475b1b3d7a2SStefano Zampini   /* Number sequential problems */
476b1b3d7a2SStefano Zampini   auxglobal_sequential = 0;
477b1b3d7a2SStefano Zampini   ierr = PCBDDCSubsetNumbering(PetscObjectComm((PetscObject)graph->l2gmap),graph->l2gmap,n_local_sequential_problems,auxlocal_sequential,PETSC_NULL,&n_sequential_problems,&auxglobal_sequential);CHKERRQ(ierr);
478b1b3d7a2SStefano Zampini 
479b1b3d7a2SStefano Zampini   /* update info in sub_schurs */
480*883469d8SStefano Zampini   if (sub_schurs->use_mumps && A) {
4811e9c79c2SStefano Zampini     ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
4821e9c79c2SStefano Zampini     sub_schurs->A = A;
4831e9c79c2SStefano Zampini   }
484b1b3d7a2SStefano Zampini   ierr = PetscObjectReference((PetscObject)S);CHKERRQ(ierr);
485b1b3d7a2SStefano Zampini   sub_schurs->S = S;
486b1b3d7a2SStefano Zampini   ierr = PetscObjectReference((PetscObject)is_I);CHKERRQ(ierr);
487b1b3d7a2SStefano Zampini   sub_schurs->is_I = is_I;
488b1b3d7a2SStefano Zampini   ierr = PetscObjectReference((PetscObject)is_B);CHKERRQ(ierr);
489b1b3d7a2SStefano Zampini   sub_schurs->is_B = is_B;
4905db18549SStefano Zampini   ierr = PetscObjectReference((PetscObject)graph->l2gmap);CHKERRQ(ierr);
4915db18549SStefano Zampini   sub_schurs->l2gmap = graph->l2gmap;
4925db18549SStefano Zampini   ierr = PetscObjectReference((PetscObject)BtoNmap);CHKERRQ(ierr);
4935db18549SStefano Zampini   sub_schurs->BtoNmap = BtoNmap;
494b1b3d7a2SStefano Zampini   sub_schurs->n_subs_seq = n_local_sequential_problems;
495b1b3d7a2SStefano Zampini   sub_schurs->n_subs_par = n_local_parallel_problems;
496b1b3d7a2SStefano Zampini   sub_schurs->n_subs_seq_g = n_sequential_problems;
497b1b3d7a2SStefano Zampini   sub_schurs->n_subs_par_g = n_parallel_problems;
498b1b3d7a2SStefano Zampini   sub_schurs->n_subs = sub_schurs->n_subs_seq + sub_schurs->n_subs_par;
499b1b3d7a2SStefano Zampini   sub_schurs->is_subs = all_cc;
500b1b3d7a2SStefano Zampini   sub_schurs->index_sequential = index_sequential;
501b1b3d7a2SStefano Zampini   sub_schurs->index_parallel = index_parallel;
502b1b3d7a2SStefano Zampini   sub_schurs->auxglobal_sequential = auxglobal_sequential;
503b1b3d7a2SStefano Zampini   sub_schurs->auxglobal_parallel = auxglobal_parallel;
504b1b3d7a2SStefano Zampini 
505b1b3d7a2SStefano Zampini   /* free workspace */
506b1b3d7a2SStefano Zampini   ierr = PetscFree(auxmapping);CHKERRQ(ierr);
507b1b3d7a2SStefano Zampini   ierr = PetscFree2(auxlocal_sequential,auxlocal_parallel);CHKERRQ(ierr);
508b1b3d7a2SStefano Zampini 
509b1b3d7a2SStefano Zampini   PetscFunctionReturn(0);
510b1b3d7a2SStefano Zampini }
511b1b3d7a2SStefano Zampini 
512b1b3d7a2SStefano Zampini #undef __FUNCT__
51334a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCSubSchursCreate"
51434a97f8cSStefano Zampini PetscErrorCode PCBDDCSubSchursCreate(PCBDDCSubSchurs *sub_schurs)
51534a97f8cSStefano Zampini {
51634a97f8cSStefano Zampini   PCBDDCSubSchurs schurs_ctx;
51734a97f8cSStefano Zampini   PetscErrorCode  ierr;
51834a97f8cSStefano Zampini 
51934a97f8cSStefano Zampini   PetscFunctionBegin;
52034a97f8cSStefano Zampini   ierr = PetscNew(&schurs_ctx);CHKERRQ(ierr);
5215ff63025SStefano Zampini   schurs_ctx->n_subs = 0;
52234a97f8cSStefano Zampini   *sub_schurs = schurs_ctx;
52334a97f8cSStefano Zampini   PetscFunctionReturn(0);
52434a97f8cSStefano Zampini }
52534a97f8cSStefano Zampini 
52634a97f8cSStefano Zampini #undef __FUNCT__
52734a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCSubSchursDestroy"
52834a97f8cSStefano Zampini PetscErrorCode PCBDDCSubSchursDestroy(PCBDDCSubSchurs *sub_schurs)
52934a97f8cSStefano Zampini {
53034a97f8cSStefano Zampini   PetscErrorCode ierr;
53134a97f8cSStefano Zampini 
53234a97f8cSStefano Zampini   PetscFunctionBegin;
53334a97f8cSStefano Zampini   ierr = PCBDDCSubSchursReset(*sub_schurs);CHKERRQ(ierr);
53434a97f8cSStefano Zampini   ierr = PetscFree(*sub_schurs);CHKERRQ(ierr);
53534a97f8cSStefano Zampini   PetscFunctionReturn(0);
53634a97f8cSStefano Zampini }
53734a97f8cSStefano Zampini 
53834a97f8cSStefano Zampini #undef __FUNCT__
53934a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCSubSchursReset"
54034a97f8cSStefano Zampini PetscErrorCode PCBDDCSubSchursReset(PCBDDCSubSchurs sub_schurs)
54134a97f8cSStefano Zampini {
54234a97f8cSStefano Zampini   PetscInt       i;
54334a97f8cSStefano Zampini   PetscErrorCode ierr;
54434a97f8cSStefano Zampini 
54534a97f8cSStefano Zampini   PetscFunctionBegin;
5461e9c79c2SStefano Zampini   ierr = MatDestroy(&sub_schurs->A);CHKERRQ(ierr);
547b1b3d7a2SStefano Zampini   ierr = MatDestroy(&sub_schurs->S);CHKERRQ(ierr);
548b1b3d7a2SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_I);CHKERRQ(ierr);
549b1b3d7a2SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_B);CHKERRQ(ierr);
5505db18549SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->l2gmap);CHKERRQ(ierr);
5515db18549SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->BtoNmap);CHKERRQ(ierr);
55241c3ba1bSStefano Zampini   ierr = MatDestroy(&sub_schurs->S_Ej_all);CHKERRQ(ierr);
55341c3ba1bSStefano Zampini   ierr = MatDestroy(&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
5545db18549SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_Ej_all);CHKERRQ(ierr);
55534a97f8cSStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
556b1b3d7a2SStefano Zampini     ierr = ISDestroy(&sub_schurs->is_subs[i]);CHKERRQ(ierr);
55734a97f8cSStefano Zampini     ierr = ISDestroy(&sub_schurs->is_AEj_B[i]);CHKERRQ(ierr);
55834a97f8cSStefano Zampini     ierr = MatDestroy(&sub_schurs->S_Ej[i]);CHKERRQ(ierr);
5599221af80SStefano Zampini     ierr = VecDestroy(&sub_schurs->work1[i]);CHKERRQ(ierr);
5609221af80SStefano Zampini     ierr = VecDestroy(&sub_schurs->work2[i]);CHKERRQ(ierr);
56134a97f8cSStefano Zampini   }
56268270318SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_I_layer);CHKERRQ(ierr);
5635ff63025SStefano Zampini   if (sub_schurs->n_subs) {
564b1b3d7a2SStefano Zampini     ierr = PetscFree(sub_schurs->is_subs);CHKERRQ(ierr);
56568270318SStefano Zampini     ierr = PetscFree4(sub_schurs->is_AEj_B,sub_schurs->S_Ej,sub_schurs->work1,sub_schurs->work2);CHKERRQ(ierr);
566b1b3d7a2SStefano Zampini     ierr = PetscFree2(sub_schurs->index_sequential,sub_schurs->index_parallel);CHKERRQ(ierr);
567b1b3d7a2SStefano Zampini     ierr = PetscFree(sub_schurs->auxglobal_sequential);CHKERRQ(ierr);
568b1b3d7a2SStefano Zampini     ierr = PetscFree(sub_schurs->auxglobal_parallel);CHKERRQ(ierr);
5695ff63025SStefano Zampini   }
57034a97f8cSStefano Zampini   sub_schurs->n_subs = 0;
57134a97f8cSStefano Zampini   PetscFunctionReturn(0);
57234a97f8cSStefano Zampini }
57334a97f8cSStefano Zampini 
57434a97f8cSStefano Zampini #undef __FUNCT__
57534a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCAdjGetNextLayer_Private"
5762a155e38SStefano Zampini PETSC_STATIC_INLINE PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt* queue_tip,PetscInt n_prev,PetscBT touched,PetscInt* xadj,PetscInt* adjncy,PetscInt* n_added)
57734a97f8cSStefano Zampini {
57834a97f8cSStefano Zampini   PetscInt       i,j,n;
57934a97f8cSStefano Zampini   PetscErrorCode ierr;
58034a97f8cSStefano Zampini 
58134a97f8cSStefano Zampini   PetscFunctionBegin;
58234a97f8cSStefano Zampini   n = 0;
58334a97f8cSStefano Zampini   for (i=-n_prev;i<0;i++) {
58434a97f8cSStefano Zampini     PetscInt start_dof = queue_tip[i];
58534a97f8cSStefano Zampini     for (j=xadj[start_dof];j<xadj[start_dof+1];j++) {
58634a97f8cSStefano Zampini       PetscInt dof = adjncy[j];
58734a97f8cSStefano Zampini       if (!PetscBTLookup(touched,dof)) {
58834a97f8cSStefano Zampini         ierr = PetscBTSet(touched,dof);CHKERRQ(ierr);
58934a97f8cSStefano Zampini         queue_tip[n] = dof;
59034a97f8cSStefano Zampini         n++;
59134a97f8cSStefano Zampini       }
59234a97f8cSStefano Zampini     }
59334a97f8cSStefano Zampini   }
59434a97f8cSStefano Zampini   *n_added = n;
59534a97f8cSStefano Zampini   PetscFunctionReturn(0);
59634a97f8cSStefano Zampini }
59734a97f8cSStefano Zampini 
59834a97f8cSStefano Zampini #undef __FUNCT__
59934a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCSubSchursSetUp"
6007821e9e7SStefano Zampini PetscErrorCode PCBDDCSubSchursSetUp(PCBDDCSubSchurs sub_schurs, Mat S, IS is_A_I, IS is_A_B, PetscInt ncc, IS is_cc[], PetscInt xadj[], PetscInt adjncy[], PetscInt nlayers)
60134a97f8cSStefano Zampini {
60234a97f8cSStefano Zampini   Mat                    A_II,A_IB,A_BI,A_BB;
6032a155e38SStefano Zampini   Mat                    AE_II,*AE_IE,*AE_EI,*AE_EE;
6042a155e38SStefano Zampini   IS                     is_I,*is_subset_B;
6052a155e38SStefano Zampini   ISLocalToGlobalMapping BtoNmap;
6062a155e38SStefano Zampini   PetscInt               i;
60734a97f8cSStefano Zampini   PetscBool              is_sorted;
60834a97f8cSStefano Zampini   PetscErrorCode         ierr;
60934a97f8cSStefano Zampini 
61034a97f8cSStefano Zampini   PetscFunctionBegin;
61134a97f8cSStefano Zampini   ierr = ISSorted(is_A_I,&is_sorted);CHKERRQ(ierr);
61234a97f8cSStefano Zampini   if (!is_sorted) {
61334a97f8cSStefano Zampini     SETERRQ(PetscObjectComm((PetscObject)is_A_I),PETSC_ERR_PLIB,"IS for I dofs should be shorted");
61434a97f8cSStefano Zampini   }
61534a97f8cSStefano Zampini   ierr = ISSorted(is_A_B,&is_sorted);CHKERRQ(ierr);
61634a97f8cSStefano Zampini   if (!is_sorted) {
61734a97f8cSStefano Zampini     SETERRQ(PetscObjectComm((PetscObject)is_A_B),PETSC_ERR_PLIB,"IS for B dofs should be shorted");
61834a97f8cSStefano Zampini   }
61934a97f8cSStefano Zampini 
62034a97f8cSStefano Zampini   /* get Schur complement matrices */
62134a97f8cSStefano Zampini   ierr = MatSchurComplementGetSubMatrices(S,&A_II,NULL,&A_IB,&A_BI,&A_BB);CHKERRQ(ierr);
62234a97f8cSStefano Zampini 
62334a97f8cSStefano Zampini   /* allocate space for schur complements */
62468270318SStefano Zampini   sub_schurs->n_subs = ncc;
62568270318SStefano Zampini   ierr = PetscMalloc4(sub_schurs->n_subs,&sub_schurs->is_AEj_B,
62641c3ba1bSStefano Zampini                       sub_schurs->n_subs,&sub_schurs->S_Ej,
62741c3ba1bSStefano Zampini                       sub_schurs->n_subs,&sub_schurs->work1,
62841c3ba1bSStefano Zampini                       sub_schurs->n_subs,&sub_schurs->work2);CHKERRQ(ierr);
6292a155e38SStefano Zampini   ierr = PetscMalloc4(ncc,&is_subset_B,ncc,&AE_IE,ncc,&AE_EI,ncc,&AE_EE);CHKERRQ(ierr);
63034a97f8cSStefano Zampini 
6312a155e38SStefano Zampini   /* maps */
6322a155e38SStefano Zampini   if (sub_schurs->n_subs && nlayers >= 0 && xadj != NULL && adjncy != NULL) { /* Interior problems can be different from the original one */
6332a155e38SStefano Zampini     ISLocalToGlobalMapping ItoNmap;
6342a155e38SStefano Zampini     PetscBT                touched;
63534a97f8cSStefano Zampini     const PetscInt*        idx_B;
6362a155e38SStefano Zampini     PetscInt               n_I,n_B,n_local_dofs,n_prev_added,j,layer,*local_numbering;
6372a155e38SStefano Zampini 
6382a155e38SStefano Zampini     /* get sizes */
6392a155e38SStefano Zampini     ierr = ISGetLocalSize(is_A_I,&n_I);CHKERRQ(ierr);
6402a155e38SStefano Zampini     ierr = ISGetLocalSize(is_A_B,&n_B);CHKERRQ(ierr);
6412a155e38SStefano Zampini 
6422a155e38SStefano Zampini     ierr = ISLocalToGlobalMappingCreateIS(is_A_I,&ItoNmap);CHKERRQ(ierr);
6432a155e38SStefano Zampini     ierr = PetscMalloc1(n_I+n_B,&local_numbering);CHKERRQ(ierr);
6442a155e38SStefano Zampini     ierr = PetscBTCreate(n_I+n_B,&touched);CHKERRQ(ierr);
6452a155e38SStefano Zampini     ierr = PetscBTMemzero(n_I+n_B,touched);CHKERRQ(ierr);
64634a97f8cSStefano Zampini 
64734a97f8cSStefano Zampini     /* all boundary dofs must be skipped when adding layers */
64834a97f8cSStefano Zampini     ierr = ISGetIndices(is_A_B,&idx_B);CHKERRQ(ierr);
64934a97f8cSStefano Zampini     for (j=0;j<n_B;j++) {
65034a97f8cSStefano Zampini       ierr = PetscBTSet(touched,idx_B[j]);CHKERRQ(ierr);
65134a97f8cSStefano Zampini     }
6522a155e38SStefano Zampini     ierr = PetscMemcpy(local_numbering,idx_B,n_B*sizeof(PetscInt));CHKERRQ(ierr);
65334a97f8cSStefano Zampini     ierr = ISRestoreIndices(is_A_B,&idx_B);CHKERRQ(ierr);
65434a97f8cSStefano Zampini 
65534a97f8cSStefano Zampini     /* add next layers of dofs */
6562a155e38SStefano Zampini     n_local_dofs = n_B;
6572a155e38SStefano Zampini     n_prev_added = n_B;
65834a97f8cSStefano Zampini     for (layer=0;layer<nlayers;layer++) {
65934a97f8cSStefano Zampini       PetscInt n_added;
6602a155e38SStefano Zampini       if (n_local_dofs == n_I+n_B) break;
6612a155e38SStefano Zampini       if (n_local_dofs > n_I+n_B) {
6622a155e38SStefano 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);
66334a97f8cSStefano Zampini       }
66434a97f8cSStefano Zampini       ierr = PCBDDCAdjGetNextLayer_Private(local_numbering+n_local_dofs,n_prev_added,touched,xadj,adjncy,&n_added);CHKERRQ(ierr);
66534a97f8cSStefano Zampini       n_prev_added = n_added;
66634a97f8cSStefano Zampini       n_local_dofs += n_added;
6679b0e569cSStefano Zampini       if (!n_added) break;
66834a97f8cSStefano Zampini     }
6692a155e38SStefano Zampini     ierr = PetscBTDestroy(&touched);CHKERRQ(ierr);
67034a97f8cSStefano Zampini 
67134a97f8cSStefano Zampini     /* IS for I dofs in original numbering and in I numbering */
67268270318SStefano Zampini     ierr = ISCreateGeneral(PetscObjectComm((PetscObject)ItoNmap),n_local_dofs-n_B,local_numbering+n_B,PETSC_COPY_VALUES,&sub_schurs->is_I_layer);CHKERRQ(ierr);
6732a155e38SStefano Zampini     ierr = PetscFree(local_numbering);CHKERRQ(ierr);
67468270318SStefano Zampini     ierr = ISSort(sub_schurs->is_I_layer);CHKERRQ(ierr);
67568270318SStefano Zampini     ierr = ISGlobalToLocalMappingApplyIS(ItoNmap,IS_GTOLM_DROP,sub_schurs->is_I_layer,&is_I);CHKERRQ(ierr);
6762a155e38SStefano Zampini     ierr = ISLocalToGlobalMappingDestroy(&ItoNmap);CHKERRQ(ierr);
67734a97f8cSStefano Zampini 
67834a97f8cSStefano Zampini     /* II block */
67934a97f8cSStefano Zampini     ierr = MatGetSubMatrix(A_II,is_I,is_I,MAT_INITIAL_MATRIX,&AE_II);CHKERRQ(ierr);
6802a155e38SStefano Zampini   } else {
6812a155e38SStefano Zampini     PetscInt n_I;
6822a155e38SStefano Zampini 
68334a97f8cSStefano Zampini     /* IS for I dofs in original numbering */
68434a97f8cSStefano Zampini     ierr = PetscObjectReference((PetscObject)is_A_I);CHKERRQ(ierr);
68568270318SStefano Zampini     sub_schurs->is_I_layer = is_A_I;
68634a97f8cSStefano Zampini 
6872a155e38SStefano Zampini     /* IS for I dofs in I numbering (strided 1) */
6882a155e38SStefano Zampini     ierr = ISGetSize(is_A_I,&n_I);CHKERRQ(ierr);
68934a97f8cSStefano Zampini     ierr = ISCreateStride(PetscObjectComm((PetscObject)is_A_I),n_I,0,1,&is_I);CHKERRQ(ierr);
69034a97f8cSStefano Zampini 
69134a97f8cSStefano Zampini     /* II block is the same */
69234a97f8cSStefano Zampini     ierr = PetscObjectReference((PetscObject)A_II);CHKERRQ(ierr);
69334a97f8cSStefano Zampini     AE_II = A_II;
69434a97f8cSStefano Zampini   }
69534a97f8cSStefano Zampini 
6962a155e38SStefano Zampini   /* subsets in original and boundary numbering */
6972a155e38SStefano Zampini   ierr = ISLocalToGlobalMappingCreateIS(is_A_B,&BtoNmap);CHKERRQ(ierr);
6982a155e38SStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
6992a155e38SStefano Zampini     ierr = ISDuplicate(is_cc[i],&sub_schurs->is_AEj_B[i]);CHKERRQ(ierr);
7002a155e38SStefano Zampini     ierr = ISSort(sub_schurs->is_AEj_B[i]);CHKERRQ(ierr);
7012a155e38SStefano Zampini     ierr = ISGlobalToLocalMappingApplyIS(BtoNmap,IS_GTOLM_DROP,sub_schurs->is_AEj_B[i],&is_subset_B[i]);CHKERRQ(ierr);
7022a155e38SStefano Zampini   }
7032a155e38SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&BtoNmap);CHKERRQ(ierr);
7042a155e38SStefano Zampini 
7052a155e38SStefano Zampini   /* EE block */
7062a155e38SStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
7072a155e38SStefano Zampini     ierr = MatGetSubMatrix(A_BB,is_subset_B[i],is_subset_B[i],MAT_INITIAL_MATRIX,&AE_EE[i]);CHKERRQ(ierr);
7082a155e38SStefano Zampini   }
7092a155e38SStefano Zampini   /* IE block */
7102a155e38SStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
7112a155e38SStefano Zampini     ierr = MatGetSubMatrix(A_IB,is_I,is_subset_B[i],MAT_INITIAL_MATRIX,&AE_IE[i]);CHKERRQ(ierr);
7122a155e38SStefano Zampini   }
71334a97f8cSStefano Zampini   /* EI block */
7142a155e38SStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
7152a155e38SStefano Zampini     ierr = MatGetSubMatrix(A_BI,is_subset_B[i],is_I,MAT_INITIAL_MATRIX,&AE_EI[i]);CHKERRQ(ierr);
7162a155e38SStefano Zampini   }
71734a97f8cSStefano Zampini 
71834a97f8cSStefano Zampini   /* setup Schur complements on subset */
7192a155e38SStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
7202a155e38SStefano Zampini     ierr = MatCreateSchurComplement(AE_II,AE_II,AE_IE[i],AE_EI[i],AE_EE[i],&sub_schurs->S_Ej[i]);CHKERRQ(ierr);
7218a26ef87SStefano Zampini     ierr = MatCreateVecs(sub_schurs->S_Ej[i],&sub_schurs->work1[i],&sub_schurs->work2[i]);CHKERRQ(ierr);
72234a97f8cSStefano Zampini     if (AE_II == A_II) { /* we can reuse the same ksp */
72334a97f8cSStefano Zampini       KSP ksp;
72434a97f8cSStefano Zampini       ierr = MatSchurComplementGetKSP(S,&ksp);CHKERRQ(ierr);
72534a97f8cSStefano Zampini       ierr = MatSchurComplementSetKSP(sub_schurs->S_Ej[i],ksp);CHKERRQ(ierr);
72634a97f8cSStefano Zampini     } else { /* build new ksp object which inherits ksp and pc types from the original one */
72734a97f8cSStefano Zampini       KSP      origksp,schurksp;
72834a97f8cSStefano Zampini       PC       origpc,schurpc;
72934a97f8cSStefano Zampini       KSPType  ksp_type;
73034a97f8cSStefano Zampini       PCType   pc_type;
73134a97f8cSStefano Zampini       PetscInt n_internal;
73234a97f8cSStefano Zampini 
73334a97f8cSStefano Zampini       ierr = MatSchurComplementGetKSP(S,&origksp);CHKERRQ(ierr);
73434a97f8cSStefano Zampini       ierr = MatSchurComplementGetKSP(sub_schurs->S_Ej[i],&schurksp);CHKERRQ(ierr);
73534a97f8cSStefano Zampini       ierr = KSPGetType(origksp,&ksp_type);CHKERRQ(ierr);
73634a97f8cSStefano Zampini       ierr = KSPSetType(schurksp,ksp_type);CHKERRQ(ierr);
73734a97f8cSStefano Zampini       ierr = KSPGetPC(schurksp,&schurpc);CHKERRQ(ierr);
73834a97f8cSStefano Zampini       ierr = KSPGetPC(origksp,&origpc);CHKERRQ(ierr);
73934a97f8cSStefano Zampini       ierr = PCGetType(origpc,&pc_type);CHKERRQ(ierr);
74034a97f8cSStefano Zampini       ierr = PCSetType(schurpc,pc_type);CHKERRQ(ierr);
74134a97f8cSStefano Zampini       ierr = ISGetSize(is_I,&n_internal);CHKERRQ(ierr);
74234a97f8cSStefano Zampini       if (n_internal) { /* UMFPACK gives error with 0 sized problems */
74334a97f8cSStefano Zampini         MatSolverPackage solver=NULL;
74434a97f8cSStefano Zampini         ierr = PCFactorGetMatSolverPackage(origpc,(const MatSolverPackage*)&solver);CHKERRQ(ierr);
74534a97f8cSStefano Zampini         if (solver) {
74634a97f8cSStefano Zampini           ierr = PCFactorSetMatSolverPackage(schurpc,solver);CHKERRQ(ierr);
74734a97f8cSStefano Zampini         }
74834a97f8cSStefano Zampini       }
74934a97f8cSStefano Zampini       ierr = KSPSetUp(schurksp);CHKERRQ(ierr);
75034a97f8cSStefano Zampini     }
75134a97f8cSStefano Zampini   }
75234a97f8cSStefano Zampini   /* free */
7532a155e38SStefano Zampini   ierr = ISDestroy(&is_I);CHKERRQ(ierr);
7542a155e38SStefano Zampini   ierr = MatDestroy(&AE_II);CHKERRQ(ierr);
7552a155e38SStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
7562a155e38SStefano Zampini     ierr = MatDestroy(&AE_EE[i]);CHKERRQ(ierr);
7572a155e38SStefano Zampini     ierr = MatDestroy(&AE_IE[i]);CHKERRQ(ierr);
7582a155e38SStefano Zampini     ierr = MatDestroy(&AE_EI[i]);CHKERRQ(ierr);
7592a155e38SStefano Zampini     ierr = ISDestroy(&is_subset_B[i]);CHKERRQ(ierr);
7602a155e38SStefano Zampini   }
7612a155e38SStefano Zampini   ierr = PetscFree4(is_subset_B,AE_IE,AE_EI,AE_EE);CHKERRQ(ierr);
76234a97f8cSStefano Zampini   PetscFunctionReturn(0);
76334a97f8cSStefano Zampini }
764