xref: /petsc/src/ksp/pc/impls/bddc/bddcschurs.c (revision 5db18549ef4e876dcd315358182a0a733898e86f)
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*5db18549SStefano Zampini   Mat                    global_schur_subsets,*submat_global_schur_subsets,work_mat;
13*5db18549SStefano Zampini   ISLocalToGlobalMapping l2gmap_subsets;
14*5db18549SStefano Zampini   IS                     is_I,*is_subset_B,temp_is;
15*5db18549SStefano Zampini   PetscInt               *nnz,*all_local_idx_G,*all_local_idx_B,*all_local_idx_N,*all_permutation_G;
16*5db18549SStefano Zampini   PetscInt               i,subset_size,max_subset_size;
17*5db18549SStefano Zampini   PetscInt               local_size,global_size;
18b1b3d7a2SStefano Zampini   PetscBool              implicit_schurs;
19b1b3d7a2SStefano Zampini   PetscErrorCode         ierr;
20b1b3d7a2SStefano Zampini 
21b1b3d7a2SStefano Zampini   PetscFunctionBegin;
22*5db18549SStefano Zampini   //ierr = PetscObjectTypeCompare((PetscObject)sub_schurs->S,MATSCHURCOMPLEMENT,&implicit_schurs);CHKERRQ(ierr);
23*5db18549SStefano Zampini   implicit_schurs = PETSC_TRUE;
24b1b3d7a2SStefano Zampini   /* allocate space for schur complements */
25b1b3d7a2SStefano Zampini   ierr = PetscMalloc5(sub_schurs->n_subs,&sub_schurs->is_AEj_I,
26b1b3d7a2SStefano Zampini                       sub_schurs->n_subs,&sub_schurs->is_AEj_B,
27b1b3d7a2SStefano Zampini                       sub_schurs->n_subs,&sub_schurs->S_Ej,
28b1b3d7a2SStefano Zampini                       sub_schurs->n_subs,&sub_schurs->work1,
29b1b3d7a2SStefano Zampini                       sub_schurs->n_subs,&sub_schurs->work2);CHKERRQ(ierr);
30b1b3d7a2SStefano Zampini 
31b1b3d7a2SStefano Zampini   /* get Schur complement matrices */
32b1b3d7a2SStefano Zampini   if (implicit_schurs) {
33b1b3d7a2SStefano Zampini     ierr = MatSchurComplementGetSubMatrices(sub_schurs->S,&A_II,NULL,&A_IB,&A_BI,&A_BB);CHKERRQ(ierr);
34b1b3d7a2SStefano Zampini     ierr = PetscMalloc4(sub_schurs->n_subs,&is_subset_B,
35b1b3d7a2SStefano Zampini                         sub_schurs->n_subs,&AE_IE,
36b1b3d7a2SStefano Zampini                         sub_schurs->n_subs,&AE_EI,
37b1b3d7a2SStefano Zampini                         sub_schurs->n_subs,&AE_EE);CHKERRQ(ierr);
38b1b3d7a2SStefano Zampini   }
39b1b3d7a2SStefano Zampini 
40b1b3d7a2SStefano Zampini   /* determine interior problems */
41b1b3d7a2SStefano Zampini   if (nlayers >= 0 && xadj != NULL && adjncy != NULL) { /* Interior problems can be different from the original one */
42b1b3d7a2SStefano Zampini     PetscBT                touched;
43b1b3d7a2SStefano Zampini     const PetscInt*        idx_B;
44b1b3d7a2SStefano Zampini     PetscInt               n_I,n_B,n_local_dofs,n_prev_added,j,layer,*local_numbering;
45b1b3d7a2SStefano Zampini 
46b1b3d7a2SStefano Zampini     /* get sizes */
47b1b3d7a2SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr);
48b1b3d7a2SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_B,&n_B);CHKERRQ(ierr);
49b1b3d7a2SStefano Zampini 
50b1b3d7a2SStefano Zampini     ierr = PetscMalloc1(n_I+n_B,&local_numbering);CHKERRQ(ierr);
51b1b3d7a2SStefano Zampini     ierr = PetscBTCreate(n_I+n_B,&touched);CHKERRQ(ierr);
52b1b3d7a2SStefano Zampini     ierr = PetscBTMemzero(n_I+n_B,touched);CHKERRQ(ierr);
53b1b3d7a2SStefano Zampini 
54b1b3d7a2SStefano Zampini     /* all boundary dofs must be skipped when adding layers */
55b1b3d7a2SStefano Zampini     ierr = ISGetIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr);
56b1b3d7a2SStefano Zampini     for (j=0;j<n_B;j++) {
57b1b3d7a2SStefano Zampini       ierr = PetscBTSet(touched,idx_B[j]);CHKERRQ(ierr);
58b1b3d7a2SStefano Zampini     }
59b1b3d7a2SStefano Zampini     ierr = PetscMemcpy(local_numbering,idx_B,n_B*sizeof(PetscInt));CHKERRQ(ierr);
60b1b3d7a2SStefano Zampini     ierr = ISRestoreIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr);
61b1b3d7a2SStefano Zampini 
62b1b3d7a2SStefano Zampini     /* add prescribed number of layers of dofs */
63b1b3d7a2SStefano Zampini     n_local_dofs = n_B;
64b1b3d7a2SStefano Zampini     n_prev_added = n_B;
65b1b3d7a2SStefano Zampini     for (layer=0;layer<nlayers;layer++) {
66b1b3d7a2SStefano Zampini       PetscInt n_added;
67b1b3d7a2SStefano Zampini       if (n_local_dofs == n_I+n_B) break;
68b1b3d7a2SStefano Zampini       if (n_local_dofs > n_I+n_B) {
69b1b3d7a2SStefano 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);
70b1b3d7a2SStefano Zampini       }
71b1b3d7a2SStefano Zampini       ierr = PCBDDCAdjGetNextLayer_Private(local_numbering+n_local_dofs,n_prev_added,touched,xadj,adjncy,&n_added);CHKERRQ(ierr);
72b1b3d7a2SStefano Zampini       n_prev_added = n_added;
73b1b3d7a2SStefano Zampini       n_local_dofs += n_added;
74b1b3d7a2SStefano Zampini       if (!n_added) break;
75b1b3d7a2SStefano Zampini     }
76b1b3d7a2SStefano Zampini     ierr = PetscBTDestroy(&touched);CHKERRQ(ierr);
77b1b3d7a2SStefano Zampini 
78b1b3d7a2SStefano Zampini     /* IS for I dofs in original numbering */
79b1b3d7a2SStefano Zampini     ierr = ISCreateGeneral(PetscObjectComm((PetscObject)sub_schurs->is_I),n_local_dofs-n_B,local_numbering+n_B,PETSC_COPY_VALUES,&sub_schurs->is_AEj_I[0]);CHKERRQ(ierr);
80b1b3d7a2SStefano Zampini     ierr = PetscFree(local_numbering);CHKERRQ(ierr);
81b1b3d7a2SStefano Zampini     ierr = ISSort(sub_schurs->is_AEj_I[0]);CHKERRQ(ierr);
82b1b3d7a2SStefano Zampini     /* IS for I dofs in boundary numbering */
83b1b3d7a2SStefano Zampini     if (implicit_schurs) {
84b1b3d7a2SStefano Zampini       ISLocalToGlobalMapping ItoNmap;
85b1b3d7a2SStefano Zampini       ierr = ISLocalToGlobalMappingCreateIS(sub_schurs->is_I,&ItoNmap);CHKERRQ(ierr);
86b1b3d7a2SStefano Zampini       ierr = ISGlobalToLocalMappingApplyIS(ItoNmap,IS_GTOLM_DROP,sub_schurs->is_AEj_I[0],&is_I);CHKERRQ(ierr);
87b1b3d7a2SStefano Zampini       ierr = ISLocalToGlobalMappingDestroy(&ItoNmap);CHKERRQ(ierr);
88b1b3d7a2SStefano Zampini 
89b1b3d7a2SStefano Zampini       /* II block */
90b1b3d7a2SStefano Zampini       ierr = MatGetSubMatrix(A_II,is_I,is_I,MAT_INITIAL_MATRIX,&AE_II);CHKERRQ(ierr);
91b1b3d7a2SStefano Zampini     }
92b1b3d7a2SStefano Zampini   } else {
93b1b3d7a2SStefano Zampini     PetscInt n_I;
94b1b3d7a2SStefano Zampini 
95b1b3d7a2SStefano Zampini     /* IS for I dofs in original numbering */
96b1b3d7a2SStefano Zampini     ierr = PetscObjectReference((PetscObject)sub_schurs->is_I);CHKERRQ(ierr);
97b1b3d7a2SStefano Zampini     sub_schurs->is_AEj_I[0] = sub_schurs->is_I;
98b1b3d7a2SStefano Zampini 
99b1b3d7a2SStefano Zampini     /* IS for I dofs in I numbering (strided 1) */
100b1b3d7a2SStefano Zampini     if (implicit_schurs) {
101b1b3d7a2SStefano Zampini       ierr = ISGetSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr);
102b1b3d7a2SStefano Zampini       ierr = ISCreateStride(PetscObjectComm((PetscObject)sub_schurs->is_I),n_I,0,1,&is_I);CHKERRQ(ierr);
103b1b3d7a2SStefano Zampini 
104b1b3d7a2SStefano Zampini       /* II block is the same */
105b1b3d7a2SStefano Zampini       ierr = PetscObjectReference((PetscObject)A_II);CHKERRQ(ierr);
106b1b3d7a2SStefano Zampini       AE_II = A_II;
107b1b3d7a2SStefano Zampini     }
108b1b3d7a2SStefano Zampini   }
109b1b3d7a2SStefano Zampini 
110b1b3d7a2SStefano Zampini   /* TODO: just for compatibility with the previous version, needs to be fixed */
111b1b3d7a2SStefano Zampini   for (i=1;i<sub_schurs->n_subs;i++) {
112b1b3d7a2SStefano Zampini     ierr = PetscObjectReference((PetscObject)sub_schurs->is_AEj_I[0]);CHKERRQ(ierr);
113b1b3d7a2SStefano Zampini     sub_schurs->is_AEj_I[i] = sub_schurs->is_AEj_I[0];
114b1b3d7a2SStefano Zampini   }
115b1b3d7a2SStefano Zampini 
116b1b3d7a2SStefano Zampini   if (implicit_schurs) {
117b1b3d7a2SStefano Zampini     /* subsets in original and boundary numbering */
118b1b3d7a2SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
119b1b3d7a2SStefano Zampini       ierr = ISDuplicate(sub_schurs->is_subs[i],&sub_schurs->is_AEj_B[i]);CHKERRQ(ierr);
120b1b3d7a2SStefano Zampini       ierr = ISSort(sub_schurs->is_AEj_B[i]);CHKERRQ(ierr);
121*5db18549SStefano Zampini       ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,sub_schurs->is_AEj_B[i],&is_subset_B[i]);CHKERRQ(ierr);
122b1b3d7a2SStefano Zampini     }
123b1b3d7a2SStefano Zampini 
124b1b3d7a2SStefano Zampini     /* EE block */
125b1b3d7a2SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
126b1b3d7a2SStefano Zampini       ierr = MatGetSubMatrix(A_BB,is_subset_B[i],is_subset_B[i],MAT_INITIAL_MATRIX,&AE_EE[i]);CHKERRQ(ierr);
127b1b3d7a2SStefano Zampini     }
128b1b3d7a2SStefano Zampini     /* IE block */
129b1b3d7a2SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
130b1b3d7a2SStefano Zampini       ierr = MatGetSubMatrix(A_IB,is_I,is_subset_B[i],MAT_INITIAL_MATRIX,&AE_IE[i]);CHKERRQ(ierr);
131b1b3d7a2SStefano Zampini     }
132b1b3d7a2SStefano Zampini     /* EI block */
133b1b3d7a2SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
134b1b3d7a2SStefano Zampini       ierr = MatGetSubMatrix(A_BI,is_subset_B[i],is_I,MAT_INITIAL_MATRIX,&AE_EI[i]);CHKERRQ(ierr);
135b1b3d7a2SStefano Zampini     }
136b1b3d7a2SStefano Zampini 
137b1b3d7a2SStefano Zampini     /* setup Schur complements on subset */
138b1b3d7a2SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
139b1b3d7a2SStefano Zampini       ierr = MatCreateSchurComplement(AE_II,AE_II,AE_IE[i],AE_EI[i],AE_EE[i],&sub_schurs->S_Ej[i]);CHKERRQ(ierr);
140b1b3d7a2SStefano Zampini       if (AE_II == A_II) { /* we can reuse the same ksp */
141b1b3d7a2SStefano Zampini         KSP ksp;
142b1b3d7a2SStefano Zampini         ierr = MatSchurComplementGetKSP(sub_schurs->S,&ksp);CHKERRQ(ierr);
143b1b3d7a2SStefano Zampini         ierr = MatSchurComplementSetKSP(sub_schurs->S_Ej[i],ksp);CHKERRQ(ierr);
144b1b3d7a2SStefano Zampini       } else { /* build new ksp object which inherits ksp and pc types from the original one */
145b1b3d7a2SStefano Zampini         KSP      origksp,schurksp;
146b1b3d7a2SStefano Zampini         PC       origpc,schurpc;
147b1b3d7a2SStefano Zampini         KSPType  ksp_type;
148b1b3d7a2SStefano Zampini         PCType   pc_type;
149b1b3d7a2SStefano Zampini         PetscInt n_internal;
150b1b3d7a2SStefano Zampini 
151b1b3d7a2SStefano Zampini         ierr = MatSchurComplementGetKSP(sub_schurs->S,&origksp);CHKERRQ(ierr);
152b1b3d7a2SStefano Zampini         ierr = MatSchurComplementGetKSP(sub_schurs->S_Ej[i],&schurksp);CHKERRQ(ierr);
153b1b3d7a2SStefano Zampini         ierr = KSPGetType(origksp,&ksp_type);CHKERRQ(ierr);
154b1b3d7a2SStefano Zampini         ierr = KSPSetType(schurksp,ksp_type);CHKERRQ(ierr);
155b1b3d7a2SStefano Zampini         ierr = KSPGetPC(schurksp,&schurpc);CHKERRQ(ierr);
156b1b3d7a2SStefano Zampini         ierr = KSPGetPC(origksp,&origpc);CHKERRQ(ierr);
157b1b3d7a2SStefano Zampini         ierr = PCGetType(origpc,&pc_type);CHKERRQ(ierr);
158b1b3d7a2SStefano Zampini         ierr = PCSetType(schurpc,pc_type);CHKERRQ(ierr);
159b1b3d7a2SStefano Zampini         ierr = ISGetSize(is_I,&n_internal);CHKERRQ(ierr);
160b1b3d7a2SStefano Zampini         if (n_internal) { /* UMFPACK gives error with 0 sized problems */
161b1b3d7a2SStefano Zampini           MatSolverPackage solver=NULL;
162b1b3d7a2SStefano Zampini           ierr = PCFactorGetMatSolverPackage(origpc,(const MatSolverPackage*)&solver);CHKERRQ(ierr);
163b1b3d7a2SStefano Zampini           if (solver) {
164b1b3d7a2SStefano Zampini             ierr = PCFactorSetMatSolverPackage(schurpc,solver);CHKERRQ(ierr);
165b1b3d7a2SStefano Zampini           }
166b1b3d7a2SStefano Zampini         }
167b1b3d7a2SStefano Zampini         ierr = KSPSetUp(schurksp);CHKERRQ(ierr);
168b1b3d7a2SStefano Zampini       }
169b1b3d7a2SStefano Zampini     }
170b1b3d7a2SStefano Zampini     /* free */
171b1b3d7a2SStefano Zampini     ierr = ISDestroy(&is_I);CHKERRQ(ierr);
172b1b3d7a2SStefano Zampini     ierr = MatDestroy(&AE_II);CHKERRQ(ierr);
173b1b3d7a2SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
174b1b3d7a2SStefano Zampini       ierr = MatDestroy(&AE_EE[i]);CHKERRQ(ierr);
175b1b3d7a2SStefano Zampini       ierr = MatDestroy(&AE_IE[i]);CHKERRQ(ierr);
176b1b3d7a2SStefano Zampini       ierr = MatDestroy(&AE_EI[i]);CHKERRQ(ierr);
177b1b3d7a2SStefano Zampini       ierr = ISDestroy(&is_subset_B[i]);CHKERRQ(ierr);
178b1b3d7a2SStefano Zampini     }
179b1b3d7a2SStefano Zampini     ierr = PetscFree4(is_subset_B,AE_IE,AE_EI,AE_EE);CHKERRQ(ierr);
180b1b3d7a2SStefano Zampini   }
181*5db18549SStefano Zampini 
182*5db18549SStefano Zampini #if 0
183*5db18549SStefano Zampini ?????
184*5db18549SStefano Zampini   if (!sub_schurs->n_subs_seq_g) {
185*5db18549SStefano Zampini     PetscFunctionReturn(0);
186*5db18549SStefano Zampini   }
187*5db18549SStefano Zampini ?????
188*5db18549SStefano Zampini #endif
189*5db18549SStefano Zampini 
190*5db18549SStefano Zampini 
191*5db18549SStefano Zampini   /* Get info on subset sizes and sum of all subsets sizes */
192*5db18549SStefano Zampini   max_subset_size = 0;
193*5db18549SStefano Zampini   local_size = 0;
194*5db18549SStefano Zampini   for (i=0;i<sub_schurs->n_subs_seq;i++) {
195*5db18549SStefano Zampini     PetscInt j = sub_schurs->index_sequential[i];
196*5db18549SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_AEj_B[j],&subset_size);CHKERRQ(ierr);
197*5db18549SStefano Zampini     max_subset_size = PetscMax(subset_size,max_subset_size);
198*5db18549SStefano Zampini     local_size += subset_size;
199*5db18549SStefano Zampini   }
200*5db18549SStefano Zampini 
201*5db18549SStefano Zampini   /* Work arrays for local indices */
202*5db18549SStefano Zampini   ierr = PetscMalloc1(local_size,&all_local_idx_B);CHKERRQ(ierr);
203*5db18549SStefano Zampini   ierr = PetscMalloc1(local_size,&all_local_idx_N);CHKERRQ(ierr);
204*5db18549SStefano Zampini   ierr = PetscMalloc1(local_size,&nnz);CHKERRQ(ierr);
205*5db18549SStefano Zampini 
206*5db18549SStefano Zampini   /* Get local indices in local whole numbering and local boundary numbering */
207*5db18549SStefano Zampini   local_size = 0;
208*5db18549SStefano Zampini   for (i=0;i<sub_schurs->n_subs_seq;i++) {
209*5db18549SStefano Zampini     PetscInt j;
210*5db18549SStefano Zampini     const    PetscInt *idxs;
211*5db18549SStefano Zampini 
212*5db18549SStefano Zampini     PetscInt local_problem_index = sub_schurs->index_sequential[i];
213*5db18549SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_AEj_B[local_problem_index],&subset_size);CHKERRQ(ierr);
214*5db18549SStefano Zampini     ierr = ISGetIndices(sub_schurs->is_AEj_B[local_problem_index],&idxs);CHKERRQ(ierr);
215*5db18549SStefano Zampini     /* subset indices in local numbering */
216*5db18549SStefano Zampini     ierr = PetscMemcpy(all_local_idx_N+local_size,idxs,subset_size*sizeof(PetscInt));CHKERRQ(ierr);
217*5db18549SStefano Zampini     ierr = ISRestoreIndices(sub_schurs->is_AEj_B[local_problem_index],&idxs);CHKERRQ(ierr);
218*5db18549SStefano Zampini     for (j=0;j<subset_size;j++) nnz[local_size+j] = subset_size;
219*5db18549SStefano Zampini     local_size += subset_size;
220*5db18549SStefano Zampini   }
221*5db18549SStefano Zampini 
222*5db18549SStefano Zampini   /* subset indices in local boundary numbering */
223*5db18549SStefano Zampini   ierr = ISGlobalToLocalMappingApply(sub_schurs->BtoNmap,IS_GTOLM_DROP,local_size,all_local_idx_N,&subset_size,all_local_idx_B);CHKERRQ(ierr);
224*5db18549SStefano Zampini   if (subset_size != local_size) {
225*5db18549SStefano Zampini     SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in sub_schurs serial (BtoNmap)! %d != %d\n",subset_size,local_size);
226*5db18549SStefano Zampini   }
227*5db18549SStefano Zampini   ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size,all_local_idx_B,PETSC_OWN_POINTER,&sub_schurs->is_Ej_all);CHKERRQ(ierr);
228*5db18549SStefano Zampini 
229*5db18549SStefano Zampini   /* Number dofs on all subsets (parallel) and sort numbering */
230*5db18549SStefano Zampini   ierr = PCBDDCSubsetNumbering(PetscObjectComm((PetscObject)sub_schurs->l2gmap),sub_schurs->l2gmap,local_size,all_local_idx_N,PETSC_NULL,&global_size,&all_local_idx_G);CHKERRQ(ierr);
231*5db18549SStefano Zampini   ierr = PetscMalloc1(local_size,&all_permutation_G);CHKERRQ(ierr);
232*5db18549SStefano Zampini   for (i=0;i<local_size;i++) {
233*5db18549SStefano Zampini     all_permutation_G[i]=i;
234*5db18549SStefano Zampini   }
235*5db18549SStefano Zampini   ierr = PetscSortIntWithPermutation(local_size,all_local_idx_G,all_permutation_G);CHKERRQ(ierr);
236*5db18549SStefano Zampini 
237*5db18549SStefano Zampini   /* Local matrix of all local Schur on subsets */
238*5db18549SStefano Zampini   ierr = MatCreate(PETSC_COMM_SELF,&sub_schurs->S_Ej_all);CHKERRQ(ierr);
239*5db18549SStefano Zampini   ierr = MatSetSizes(sub_schurs->S_Ej_all,PETSC_DECIDE,PETSC_DECIDE,local_size,local_size);CHKERRQ(ierr);
240*5db18549SStefano Zampini   ierr = MatSetType(sub_schurs->S_Ej_all,MATAIJ);CHKERRQ(ierr);
241*5db18549SStefano Zampini   ierr = MatSeqAIJSetPreallocation(sub_schurs->S_Ej_all,0,nnz);CHKERRQ(ierr);
242*5db18549SStefano Zampini   ierr = PetscFree(nnz);CHKERRQ(ierr);
243*5db18549SStefano Zampini 
244*5db18549SStefano Zampini   if (implicit_schurs) {
245*5db18549SStefano Zampini     PetscScalar *fill_vals;
246*5db18549SStefano Zampini     PetscInt    *dummy_idx;
247*5db18549SStefano Zampini 
248*5db18549SStefano Zampini     ierr = MatSetOption(sub_schurs->S_Ej_all,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr);
249*5db18549SStefano Zampini 
250*5db18549SStefano Zampini     /* Work arrays */
251*5db18549SStefano Zampini     ierr = PetscMalloc2(max_subset_size,&dummy_idx,max_subset_size*max_subset_size,&fill_vals);CHKERRQ(ierr);
252*5db18549SStefano Zampini 
253*5db18549SStefano Zampini     /* Loop on local problems to compute Schur complements explicitly (TODO; optimize)*/
254*5db18549SStefano Zampini     local_size = 0;
255*5db18549SStefano Zampini     for (i=0;i<sub_schurs->n_subs_seq;i++) {
256*5db18549SStefano Zampini       Vec work1,work2;
257*5db18549SStefano Zampini       PetscInt j,local_problem_index = sub_schurs->index_sequential[i];
258*5db18549SStefano Zampini 
259*5db18549SStefano Zampini       ierr = MatCreateVecs(sub_schurs->S_Ej[i],&work1,&work2);CHKERRQ(ierr);
260*5db18549SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_AEj_B[local_problem_index],&subset_size);CHKERRQ(ierr);
261*5db18549SStefano Zampini       /* local Schur */
262*5db18549SStefano Zampini       for (j=0;j<subset_size;j++) {
263*5db18549SStefano Zampini         ierr = VecSet(work1,0.0);CHKERRQ(ierr);
264*5db18549SStefano Zampini         ierr = VecSetValue(work1,j,1.0,INSERT_VALUES);CHKERRQ(ierr);
265*5db18549SStefano Zampini         ierr = VecPlaceArray(work2,&fill_vals[j*subset_size]);CHKERRQ(ierr);
266*5db18549SStefano Zampini         ierr = MatMult(sub_schurs->S_Ej[local_problem_index],work1,work2);CHKERRQ(ierr);
267*5db18549SStefano Zampini         ierr = VecResetArray(work2);CHKERRQ(ierr);
268*5db18549SStefano Zampini       }
269*5db18549SStefano Zampini       for (j=0;j<subset_size;j++) {
270*5db18549SStefano Zampini         dummy_idx[j]=local_size+j;
271*5db18549SStefano Zampini       }
272*5db18549SStefano Zampini       ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,fill_vals,INSERT_VALUES);CHKERRQ(ierr);
273*5db18549SStefano Zampini       local_size += subset_size;
274*5db18549SStefano Zampini       ierr = VecDestroy(&work1);CHKERRQ(ierr);
275*5db18549SStefano Zampini       ierr = VecDestroy(&work2);CHKERRQ(ierr);
276*5db18549SStefano Zampini     }
277*5db18549SStefano Zampini     ierr = PetscFree2(dummy_idx,fill_vals);CHKERRQ(ierr);
278*5db18549SStefano Zampini   }
279*5db18549SStefano Zampini   ierr = MatAssemblyBegin(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
280*5db18549SStefano Zampini   ierr = MatAssemblyEnd(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
281*5db18549SStefano Zampini 
282*5db18549SStefano Zampini   /* Global matrix of all assembled Schur on subsets */
283*5db18549SStefano Zampini   ierr = ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)sub_schurs->l2gmap),1,local_size,all_local_idx_G,PETSC_COPY_VALUES,&l2gmap_subsets);CHKERRQ(ierr);
284*5db18549SStefano Zampini   ierr = MatCreateIS(PetscObjectComm((PetscObject)sub_schurs->l2gmap),1,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size,l2gmap_subsets,&work_mat);CHKERRQ(ierr);
285*5db18549SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&l2gmap_subsets);CHKERRQ(ierr);
286*5db18549SStefano Zampini   ierr = MatISSetLocalMat(work_mat,sub_schurs->S_Ej_all);CHKERRQ(ierr);
287*5db18549SStefano Zampini   ierr = MatISGetMPIXAIJ(work_mat,MAT_INITIAL_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
288*5db18549SStefano Zampini   ierr = MatDestroy(&work_mat);CHKERRQ(ierr);
289*5db18549SStefano Zampini 
290*5db18549SStefano Zampini   /* Get local part of (\sum_j S_Ej) */
291*5db18549SStefano Zampini   for (i=0;i<local_size;i++) {
292*5db18549SStefano Zampini     all_local_idx_N[i] = all_local_idx_G[all_permutation_G[i]];
293*5db18549SStefano Zampini   }
294*5db18549SStefano Zampini   ierr = ISCreateGeneral(PetscObjectComm((PetscObject)sub_schurs->l2gmap),local_size,all_local_idx_N,PETSC_OWN_POINTER,&temp_is);CHKERRQ(ierr);
295*5db18549SStefano Zampini   ierr = MatGetSubMatrices(global_schur_subsets,1,&temp_is,&temp_is,MAT_INITIAL_MATRIX,&submat_global_schur_subsets);CHKERRQ(ierr);
296*5db18549SStefano Zampini   ierr = MatDestroy(&global_schur_subsets);CHKERRQ(ierr);
297*5db18549SStefano Zampini   ierr = ISDestroy(&temp_is);CHKERRQ(ierr);
298*5db18549SStefano Zampini   for (i=0;i<local_size;i++) {
299*5db18549SStefano Zampini     all_local_idx_G[all_permutation_G[i]] = i;
300*5db18549SStefano Zampini   }
301*5db18549SStefano Zampini   ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size,all_local_idx_G,PETSC_OWN_POINTER,&temp_is);CHKERRQ(ierr);
302*5db18549SStefano Zampini   ierr = ISSetPermutation(temp_is);CHKERRQ(ierr);
303*5db18549SStefano Zampini   ierr = MatPermute(submat_global_schur_subsets[0],temp_is,temp_is,&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
304*5db18549SStefano Zampini   ierr = MatDestroyMatrices(1,&submat_global_schur_subsets);CHKERRQ(ierr);
305*5db18549SStefano Zampini   ierr = ISDestroy(&temp_is);CHKERRQ(ierr);
306*5db18549SStefano Zampini   ierr = PetscFree(all_permutation_G);CHKERRQ(ierr);
307b1b3d7a2SStefano Zampini   PetscFunctionReturn(0);
308b1b3d7a2SStefano Zampini }
309b1b3d7a2SStefano Zampini 
310b1b3d7a2SStefano Zampini #undef __FUNCT__
311b1b3d7a2SStefano Zampini #define __FUNCT__ "PCBDDCSubSchursInit"
312*5db18549SStefano Zampini PetscErrorCode PCBDDCSubSchursInit(PCBDDCSubSchurs sub_schurs, Mat A, Mat S, IS is_I, IS is_B, PCBDDCGraph graph, ISLocalToGlobalMapping BtoNmap, PetscInt seqthreshold)
313b1b3d7a2SStefano Zampini {
314b1b3d7a2SStefano Zampini   IS                  *faces,*edges,*all_cc;
315b1b3d7a2SStefano Zampini   PetscInt            *index_sequential,*index_parallel;
316b1b3d7a2SStefano Zampini   PetscInt            *auxlocal_sequential,*auxlocal_parallel;
317b1b3d7a2SStefano Zampini   PetscInt            *auxglobal_sequential,*auxglobal_parallel;
318b1b3d7a2SStefano Zampini   PetscInt            *auxmapping;//,*idxs;
319b1b3d7a2SStefano Zampini   PetscInt            i,max_subset_size;
320b1b3d7a2SStefano Zampini   PetscInt            n_sequential_problems,n_local_sequential_problems,n_parallel_problems,n_local_parallel_problems;
321b1b3d7a2SStefano Zampini   PetscInt            n_faces,n_edges,n_all_cc;
322b1b3d7a2SStefano Zampini   PetscBool           is_sorted;
323b1b3d7a2SStefano Zampini   PetscErrorCode  ierr;
324b1b3d7a2SStefano Zampini 
325b1b3d7a2SStefano Zampini   PetscFunctionBegin;
326b1b3d7a2SStefano Zampini   ierr = ISSorted(is_I,&is_sorted);CHKERRQ(ierr);
327b1b3d7a2SStefano Zampini   if (!is_sorted) {
328b1b3d7a2SStefano Zampini     SETERRQ(PetscObjectComm((PetscObject)is_I),PETSC_ERR_PLIB,"IS for I dofs should be shorted");
329b1b3d7a2SStefano Zampini   }
330b1b3d7a2SStefano Zampini   ierr = ISSorted(is_B,&is_sorted);CHKERRQ(ierr);
331b1b3d7a2SStefano Zampini   if (!is_sorted) {
332b1b3d7a2SStefano Zampini     SETERRQ(PetscObjectComm((PetscObject)is_B),PETSC_ERR_PLIB,"IS for B dofs should be shorted");
333b1b3d7a2SStefano Zampini   }
334b1b3d7a2SStefano Zampini 
335b1b3d7a2SStefano Zampini   /* reset any previous data */
336b1b3d7a2SStefano Zampini   ierr = PCBDDCSubSchursReset(sub_schurs);CHKERRQ(ierr);
337b1b3d7a2SStefano Zampini 
338b1b3d7a2SStefano Zampini   /* get index sets for faces and edges */
339b1b3d7a2SStefano Zampini   ierr = PCBDDCGraphGetCandidatesIS(graph,&n_faces,&faces,&n_edges,&edges,NULL);CHKERRQ(ierr);
340b1b3d7a2SStefano Zampini   n_all_cc = n_faces+n_edges;
341b1b3d7a2SStefano Zampini   ierr = PetscMalloc1(n_all_cc,&all_cc);CHKERRQ(ierr);
342b1b3d7a2SStefano Zampini   for (i=0;i<n_faces;i++) {
343b1b3d7a2SStefano Zampini     all_cc[i] = faces[i];
344b1b3d7a2SStefano Zampini   }
345b1b3d7a2SStefano Zampini   for (i=0;i<n_edges;i++) {
346b1b3d7a2SStefano Zampini     all_cc[n_faces+i] = edges[i];
347b1b3d7a2SStefano Zampini   }
348b1b3d7a2SStefano Zampini   ierr = PetscFree(faces);CHKERRQ(ierr);
349b1b3d7a2SStefano Zampini   ierr = PetscFree(edges);CHKERRQ(ierr);
350b1b3d7a2SStefano Zampini 
351b1b3d7a2SStefano Zampini   /* map interface's subsets */
352b1b3d7a2SStefano Zampini   max_subset_size = 0;
353b1b3d7a2SStefano Zampini   for (i=0;i<n_all_cc;i++) {
354b1b3d7a2SStefano Zampini     PetscInt subset_size;
355b1b3d7a2SStefano Zampini     ierr = ISGetLocalSize(all_cc[i],&subset_size);CHKERRQ(ierr);
356b1b3d7a2SStefano Zampini     max_subset_size = PetscMax(max_subset_size,subset_size);
357b1b3d7a2SStefano Zampini   }
358b1b3d7a2SStefano Zampini   ierr = PetscMalloc1(max_subset_size,&auxmapping);CHKERRQ(ierr);
359b1b3d7a2SStefano Zampini   ierr = PetscMalloc2(graph->ncc,&auxlocal_sequential,
360b1b3d7a2SStefano Zampini                       graph->ncc,&auxlocal_parallel);CHKERRQ(ierr);
361b1b3d7a2SStefano Zampini   ierr = PetscMalloc2(graph->ncc,&index_sequential,
362b1b3d7a2SStefano Zampini                       graph->ncc,&index_parallel);CHKERRQ(ierr);
363b1b3d7a2SStefano Zampini 
364b1b3d7a2SStefano Zampini   /* if threshold is negative, uses all sequential problems */
365b1b3d7a2SStefano Zampini   if (seqthreshold < 0) seqthreshold = max_subset_size;
366b1b3d7a2SStefano Zampini 
367b1b3d7a2SStefano Zampini   /* determine which problem has to be solved in parallel or sequentially */
368b1b3d7a2SStefano Zampini   n_local_sequential_problems = 0;
369b1b3d7a2SStefano Zampini   n_local_parallel_problems = 0;
370b1b3d7a2SStefano Zampini   for (i=0;i<n_all_cc;i++) {
371b1b3d7a2SStefano Zampini     PetscInt       subset_size,j,min_loc = 0;
372b1b3d7a2SStefano Zampini     const PetscInt *idxs;
373b1b3d7a2SStefano Zampini 
374b1b3d7a2SStefano Zampini     ierr = ISGetLocalSize(all_cc[i],&subset_size);CHKERRQ(ierr);
375b1b3d7a2SStefano Zampini     ierr = ISGetIndices(all_cc[i],&idxs);CHKERRQ(ierr);
376b1b3d7a2SStefano Zampini     ierr = ISLocalToGlobalMappingApply(graph->l2gmap,subset_size,idxs,auxmapping);CHKERRQ(ierr);
377b1b3d7a2SStefano Zampini     for (j=1;j<subset_size;j++) {
378b1b3d7a2SStefano Zampini       if (auxmapping[j]<auxmapping[min_loc]) {
379b1b3d7a2SStefano Zampini         min_loc = j;
380b1b3d7a2SStefano Zampini       }
381b1b3d7a2SStefano Zampini     }
382b1b3d7a2SStefano Zampini     if (subset_size > seqthreshold) {
383b1b3d7a2SStefano Zampini       index_parallel[n_local_parallel_problems] = i;
384b1b3d7a2SStefano Zampini       auxlocal_parallel[n_local_parallel_problems] = idxs[min_loc];
385b1b3d7a2SStefano Zampini       n_local_parallel_problems++;
386b1b3d7a2SStefano Zampini     } else {
387b1b3d7a2SStefano Zampini       index_sequential[n_local_sequential_problems] = i;
388b1b3d7a2SStefano Zampini       auxlocal_sequential[n_local_sequential_problems] = idxs[min_loc];
389b1b3d7a2SStefano Zampini       n_local_sequential_problems++;
390b1b3d7a2SStefano Zampini     }
391b1b3d7a2SStefano Zampini     ierr = ISRestoreIndices(all_cc[i],&idxs);CHKERRQ(ierr);
392b1b3d7a2SStefano Zampini   }
393b1b3d7a2SStefano Zampini 
394b1b3d7a2SStefano Zampini   /* Number parallel problems */
395b1b3d7a2SStefano Zampini   auxglobal_parallel = 0;
396b1b3d7a2SStefano Zampini   ierr = PCBDDCSubsetNumbering(PetscObjectComm((PetscObject)graph->l2gmap),graph->l2gmap,n_local_parallel_problems,auxlocal_parallel,PETSC_NULL,&n_parallel_problems,&auxglobal_parallel);CHKERRQ(ierr);
397b1b3d7a2SStefano Zampini 
398b1b3d7a2SStefano Zampini   /* Number sequential problems */
399b1b3d7a2SStefano Zampini   auxglobal_sequential = 0;
400b1b3d7a2SStefano Zampini   ierr = PCBDDCSubsetNumbering(PetscObjectComm((PetscObject)graph->l2gmap),graph->l2gmap,n_local_sequential_problems,auxlocal_sequential,PETSC_NULL,&n_sequential_problems,&auxglobal_sequential);CHKERRQ(ierr);
401b1b3d7a2SStefano Zampini 
402b1b3d7a2SStefano Zampini   /* update info in sub_schurs */
4031e9c79c2SStefano Zampini   if (A) {
4041e9c79c2SStefano Zampini     ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
4051e9c79c2SStefano Zampini     sub_schurs->A = A;
4061e9c79c2SStefano Zampini   }
407b1b3d7a2SStefano Zampini   ierr = PetscObjectReference((PetscObject)S);CHKERRQ(ierr);
408b1b3d7a2SStefano Zampini   sub_schurs->S = S;
409b1b3d7a2SStefano Zampini   ierr = PetscObjectReference((PetscObject)is_I);CHKERRQ(ierr);
410b1b3d7a2SStefano Zampini   sub_schurs->is_I = is_I;
411b1b3d7a2SStefano Zampini   ierr = PetscObjectReference((PetscObject)is_B);CHKERRQ(ierr);
412b1b3d7a2SStefano Zampini   sub_schurs->is_B = is_B;
413*5db18549SStefano Zampini   ierr = PetscObjectReference((PetscObject)graph->l2gmap);CHKERRQ(ierr);
414*5db18549SStefano Zampini   sub_schurs->l2gmap = graph->l2gmap;
415*5db18549SStefano Zampini   ierr = PetscObjectReference((PetscObject)BtoNmap);CHKERRQ(ierr);
416*5db18549SStefano Zampini   sub_schurs->BtoNmap = BtoNmap;
417b1b3d7a2SStefano Zampini   sub_schurs->n_subs_seq = n_local_sequential_problems;
418b1b3d7a2SStefano Zampini   sub_schurs->n_subs_par = n_local_parallel_problems;
419b1b3d7a2SStefano Zampini   sub_schurs->n_subs_seq_g = n_sequential_problems;
420b1b3d7a2SStefano Zampini   sub_schurs->n_subs_par_g = n_parallel_problems;
421b1b3d7a2SStefano Zampini   sub_schurs->n_subs = sub_schurs->n_subs_seq + sub_schurs->n_subs_par;
422b1b3d7a2SStefano Zampini   sub_schurs->is_subs = all_cc;
423b1b3d7a2SStefano Zampini   sub_schurs->index_sequential = index_sequential;
424b1b3d7a2SStefano Zampini   sub_schurs->index_parallel = index_parallel;
425b1b3d7a2SStefano Zampini   sub_schurs->auxglobal_sequential = auxglobal_sequential;
426b1b3d7a2SStefano Zampini   sub_schurs->auxglobal_parallel = auxglobal_parallel;
427b1b3d7a2SStefano Zampini 
428b1b3d7a2SStefano Zampini   /* free workspace */
429b1b3d7a2SStefano Zampini   ierr = PetscFree(auxmapping);CHKERRQ(ierr);
430b1b3d7a2SStefano Zampini   ierr = PetscFree2(auxlocal_sequential,auxlocal_parallel);CHKERRQ(ierr);
431b1b3d7a2SStefano Zampini 
432b1b3d7a2SStefano Zampini   PetscFunctionReturn(0);
433b1b3d7a2SStefano Zampini }
434b1b3d7a2SStefano Zampini 
435b1b3d7a2SStefano Zampini #undef __FUNCT__
43634a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCSubSchursCreate"
43734a97f8cSStefano Zampini PetscErrorCode PCBDDCSubSchursCreate(PCBDDCSubSchurs *sub_schurs)
43834a97f8cSStefano Zampini {
43934a97f8cSStefano Zampini   PCBDDCSubSchurs schurs_ctx;
44034a97f8cSStefano Zampini   PetscErrorCode  ierr;
44134a97f8cSStefano Zampini 
44234a97f8cSStefano Zampini   PetscFunctionBegin;
44334a97f8cSStefano Zampini   ierr = PetscNew(&schurs_ctx);CHKERRQ(ierr);
4445ff63025SStefano Zampini   schurs_ctx->n_subs = 0;
44534a97f8cSStefano Zampini   *sub_schurs = schurs_ctx;
44634a97f8cSStefano Zampini   PetscFunctionReturn(0);
44734a97f8cSStefano Zampini }
44834a97f8cSStefano Zampini 
44934a97f8cSStefano Zampini #undef __FUNCT__
45034a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCSubSchursDestroy"
45134a97f8cSStefano Zampini PetscErrorCode PCBDDCSubSchursDestroy(PCBDDCSubSchurs *sub_schurs)
45234a97f8cSStefano Zampini {
45334a97f8cSStefano Zampini   PetscErrorCode ierr;
45434a97f8cSStefano Zampini 
45534a97f8cSStefano Zampini   PetscFunctionBegin;
45634a97f8cSStefano Zampini   ierr = PCBDDCSubSchursReset(*sub_schurs);CHKERRQ(ierr);
45734a97f8cSStefano Zampini   ierr = PetscFree(*sub_schurs);CHKERRQ(ierr);
45834a97f8cSStefano Zampini   PetscFunctionReturn(0);
45934a97f8cSStefano Zampini }
46034a97f8cSStefano Zampini 
46134a97f8cSStefano Zampini #undef __FUNCT__
46234a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCSubSchursReset"
46334a97f8cSStefano Zampini PetscErrorCode PCBDDCSubSchursReset(PCBDDCSubSchurs sub_schurs)
46434a97f8cSStefano Zampini {
46534a97f8cSStefano Zampini   PetscInt       i;
46634a97f8cSStefano Zampini   PetscErrorCode ierr;
46734a97f8cSStefano Zampini 
46834a97f8cSStefano Zampini   PetscFunctionBegin;
4691e9c79c2SStefano Zampini   ierr = MatDestroy(&sub_schurs->A);CHKERRQ(ierr);
470b1b3d7a2SStefano Zampini   ierr = MatDestroy(&sub_schurs->S);CHKERRQ(ierr);
471b1b3d7a2SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_I);CHKERRQ(ierr);
472b1b3d7a2SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_B);CHKERRQ(ierr);
473*5db18549SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->l2gmap);CHKERRQ(ierr);
474*5db18549SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->BtoNmap);CHKERRQ(ierr);
47541c3ba1bSStefano Zampini   ierr = MatDestroy(&sub_schurs->S_Ej_all);CHKERRQ(ierr);
47641c3ba1bSStefano Zampini   ierr = MatDestroy(&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
477*5db18549SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_Ej_all);CHKERRQ(ierr);
47834a97f8cSStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
479b1b3d7a2SStefano Zampini     ierr = ISDestroy(&sub_schurs->is_subs[i]);CHKERRQ(ierr);
48034a97f8cSStefano Zampini     ierr = ISDestroy(&sub_schurs->is_AEj_I[i]);CHKERRQ(ierr);
48134a97f8cSStefano Zampini     ierr = ISDestroy(&sub_schurs->is_AEj_B[i]);CHKERRQ(ierr);
48234a97f8cSStefano Zampini     ierr = MatDestroy(&sub_schurs->S_Ej[i]);CHKERRQ(ierr);
483*5db18549SStefano Zampini     //ierr = VecDestroy(&sub_schurs->work1[i]);CHKERRQ(ierr);
484*5db18549SStefano Zampini     //ierr = VecDestroy(&sub_schurs->work2[i]);CHKERRQ(ierr);
48534a97f8cSStefano Zampini   }
4865ff63025SStefano Zampini   if (sub_schurs->n_subs) {
487b1b3d7a2SStefano Zampini     ierr = PetscFree(sub_schurs->is_subs);CHKERRQ(ierr);
48834a97f8cSStefano Zampini     ierr = PetscFree5(sub_schurs->is_AEj_I,sub_schurs->is_AEj_B,sub_schurs->S_Ej,sub_schurs->work1,sub_schurs->work2);CHKERRQ(ierr);
489b1b3d7a2SStefano Zampini     ierr = PetscFree2(sub_schurs->index_sequential,sub_schurs->index_parallel);CHKERRQ(ierr);
490b1b3d7a2SStefano Zampini     ierr = PetscFree(sub_schurs->auxglobal_sequential);CHKERRQ(ierr);
491b1b3d7a2SStefano Zampini     ierr = PetscFree(sub_schurs->auxglobal_parallel);CHKERRQ(ierr);
4925ff63025SStefano Zampini   }
49334a97f8cSStefano Zampini   sub_schurs->n_subs = 0;
49434a97f8cSStefano Zampini   PetscFunctionReturn(0);
49534a97f8cSStefano Zampini }
49634a97f8cSStefano Zampini 
49734a97f8cSStefano Zampini #undef __FUNCT__
49834a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCAdjGetNextLayer_Private"
4992a155e38SStefano Zampini PETSC_STATIC_INLINE PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt* queue_tip,PetscInt n_prev,PetscBT touched,PetscInt* xadj,PetscInt* adjncy,PetscInt* n_added)
50034a97f8cSStefano Zampini {
50134a97f8cSStefano Zampini   PetscInt       i,j,n;
50234a97f8cSStefano Zampini   PetscErrorCode ierr;
50334a97f8cSStefano Zampini 
50434a97f8cSStefano Zampini   PetscFunctionBegin;
50534a97f8cSStefano Zampini   n = 0;
50634a97f8cSStefano Zampini   for (i=-n_prev;i<0;i++) {
50734a97f8cSStefano Zampini     PetscInt start_dof = queue_tip[i];
50834a97f8cSStefano Zampini     for (j=xadj[start_dof];j<xadj[start_dof+1];j++) {
50934a97f8cSStefano Zampini       PetscInt dof = adjncy[j];
51034a97f8cSStefano Zampini       if (!PetscBTLookup(touched,dof)) {
51134a97f8cSStefano Zampini         ierr = PetscBTSet(touched,dof);CHKERRQ(ierr);
51234a97f8cSStefano Zampini         queue_tip[n] = dof;
51334a97f8cSStefano Zampini         n++;
51434a97f8cSStefano Zampini       }
51534a97f8cSStefano Zampini     }
51634a97f8cSStefano Zampini   }
51734a97f8cSStefano Zampini   *n_added = n;
51834a97f8cSStefano Zampini   PetscFunctionReturn(0);
51934a97f8cSStefano Zampini }
52034a97f8cSStefano Zampini 
52134a97f8cSStefano Zampini #undef __FUNCT__
52234a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCSubSchursSetUp"
5237821e9e7SStefano 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)
52434a97f8cSStefano Zampini {
52534a97f8cSStefano Zampini   Mat                    A_II,A_IB,A_BI,A_BB;
5262a155e38SStefano Zampini   Mat                    AE_II,*AE_IE,*AE_EI,*AE_EE;
5272a155e38SStefano Zampini   IS                     is_I,*is_subset_B;
5282a155e38SStefano Zampini   ISLocalToGlobalMapping BtoNmap;
5292a155e38SStefano Zampini   PetscInt               i;
53034a97f8cSStefano Zampini   PetscBool              is_sorted;
53134a97f8cSStefano Zampini   PetscErrorCode         ierr;
53234a97f8cSStefano Zampini 
53334a97f8cSStefano Zampini   PetscFunctionBegin;
53434a97f8cSStefano Zampini   ierr = ISSorted(is_A_I,&is_sorted);CHKERRQ(ierr);
53534a97f8cSStefano Zampini   if (!is_sorted) {
53634a97f8cSStefano Zampini     SETERRQ(PetscObjectComm((PetscObject)is_A_I),PETSC_ERR_PLIB,"IS for I dofs should be shorted");
53734a97f8cSStefano Zampini   }
53834a97f8cSStefano Zampini   ierr = ISSorted(is_A_B,&is_sorted);CHKERRQ(ierr);
53934a97f8cSStefano Zampini   if (!is_sorted) {
54034a97f8cSStefano Zampini     SETERRQ(PetscObjectComm((PetscObject)is_A_B),PETSC_ERR_PLIB,"IS for B dofs should be shorted");
54134a97f8cSStefano Zampini   }
54234a97f8cSStefano Zampini 
54334a97f8cSStefano Zampini   /* get Schur complement matrices */
54434a97f8cSStefano Zampini   ierr = MatSchurComplementGetSubMatrices(S,&A_II,NULL,&A_IB,&A_BI,&A_BB);CHKERRQ(ierr);
54534a97f8cSStefano Zampini 
54634a97f8cSStefano Zampini   /* allocate space for schur complements */
54741c3ba1bSStefano Zampini   ierr = PetscMalloc5(sub_schurs->n_subs,&sub_schurs->is_AEj_I,
54841c3ba1bSStefano Zampini                       sub_schurs->n_subs,&sub_schurs->is_AEj_B,
54941c3ba1bSStefano Zampini                       sub_schurs->n_subs,&sub_schurs->S_Ej,
55041c3ba1bSStefano Zampini                       sub_schurs->n_subs,&sub_schurs->work1,
55141c3ba1bSStefano Zampini                       sub_schurs->n_subs,&sub_schurs->work2);CHKERRQ(ierr);
5522a155e38SStefano Zampini   ierr = PetscMalloc4(ncc,&is_subset_B,ncc,&AE_IE,ncc,&AE_EI,ncc,&AE_EE);CHKERRQ(ierr);
55334a97f8cSStefano Zampini   sub_schurs->n_subs = ncc;
55434a97f8cSStefano Zampini 
5552a155e38SStefano Zampini   /* maps */
5562a155e38SStefano Zampini   if (sub_schurs->n_subs && nlayers >= 0 && xadj != NULL && adjncy != NULL) { /* Interior problems can be different from the original one */
5572a155e38SStefano Zampini     ISLocalToGlobalMapping ItoNmap;
5582a155e38SStefano Zampini     PetscBT                touched;
55934a97f8cSStefano Zampini     const PetscInt*        idx_B;
5602a155e38SStefano Zampini     PetscInt               n_I,n_B,n_local_dofs,n_prev_added,j,layer,*local_numbering;
5612a155e38SStefano Zampini 
5622a155e38SStefano Zampini     /* get sizes */
5632a155e38SStefano Zampini     ierr = ISGetLocalSize(is_A_I,&n_I);CHKERRQ(ierr);
5642a155e38SStefano Zampini     ierr = ISGetLocalSize(is_A_B,&n_B);CHKERRQ(ierr);
5652a155e38SStefano Zampini 
5662a155e38SStefano Zampini     ierr = ISLocalToGlobalMappingCreateIS(is_A_I,&ItoNmap);CHKERRQ(ierr);
5672a155e38SStefano Zampini     ierr = PetscMalloc1(n_I+n_B,&local_numbering);CHKERRQ(ierr);
5682a155e38SStefano Zampini     ierr = PetscBTCreate(n_I+n_B,&touched);CHKERRQ(ierr);
5692a155e38SStefano Zampini     ierr = PetscBTMemzero(n_I+n_B,touched);CHKERRQ(ierr);
57034a97f8cSStefano Zampini 
57134a97f8cSStefano Zampini     /* all boundary dofs must be skipped when adding layers */
57234a97f8cSStefano Zampini     ierr = ISGetIndices(is_A_B,&idx_B);CHKERRQ(ierr);
57334a97f8cSStefano Zampini     for (j=0;j<n_B;j++) {
57434a97f8cSStefano Zampini       ierr = PetscBTSet(touched,idx_B[j]);CHKERRQ(ierr);
57534a97f8cSStefano Zampini     }
5762a155e38SStefano Zampini     ierr = PetscMemcpy(local_numbering,idx_B,n_B*sizeof(PetscInt));CHKERRQ(ierr);
57734a97f8cSStefano Zampini     ierr = ISRestoreIndices(is_A_B,&idx_B);CHKERRQ(ierr);
57834a97f8cSStefano Zampini 
57934a97f8cSStefano Zampini     /* add next layers of dofs */
5802a155e38SStefano Zampini     n_local_dofs = n_B;
5812a155e38SStefano Zampini     n_prev_added = n_B;
58234a97f8cSStefano Zampini     for (layer=0;layer<nlayers;layer++) {
58334a97f8cSStefano Zampini       PetscInt n_added;
5842a155e38SStefano Zampini       if (n_local_dofs == n_I+n_B) break;
5852a155e38SStefano Zampini       if (n_local_dofs > n_I+n_B) {
5862a155e38SStefano 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);
58734a97f8cSStefano Zampini       }
58834a97f8cSStefano Zampini       ierr = PCBDDCAdjGetNextLayer_Private(local_numbering+n_local_dofs,n_prev_added,touched,xadj,adjncy,&n_added);CHKERRQ(ierr);
58934a97f8cSStefano Zampini       n_prev_added = n_added;
59034a97f8cSStefano Zampini       n_local_dofs += n_added;
5919b0e569cSStefano Zampini       if (!n_added) break;
59234a97f8cSStefano Zampini     }
5932a155e38SStefano Zampini     ierr = PetscBTDestroy(&touched);CHKERRQ(ierr);
59434a97f8cSStefano Zampini 
59534a97f8cSStefano Zampini     /* IS for I dofs in original numbering and in I numbering */
5962a155e38SStefano Zampini     ierr = ISCreateGeneral(PetscObjectComm((PetscObject)ItoNmap),n_local_dofs-n_B,local_numbering+n_B,PETSC_COPY_VALUES,&sub_schurs->is_AEj_I[0]);CHKERRQ(ierr);
5972a155e38SStefano Zampini     ierr = PetscFree(local_numbering);CHKERRQ(ierr);
5982a155e38SStefano Zampini     ierr = ISSort(sub_schurs->is_AEj_I[0]);CHKERRQ(ierr);
5992a155e38SStefano Zampini     ierr = ISGlobalToLocalMappingApplyIS(ItoNmap,IS_GTOLM_DROP,sub_schurs->is_AEj_I[0],&is_I);CHKERRQ(ierr);
6002a155e38SStefano Zampini     ierr = ISLocalToGlobalMappingDestroy(&ItoNmap);CHKERRQ(ierr);
60134a97f8cSStefano Zampini 
60234a97f8cSStefano Zampini     /* II block */
60334a97f8cSStefano Zampini     ierr = MatGetSubMatrix(A_II,is_I,is_I,MAT_INITIAL_MATRIX,&AE_II);CHKERRQ(ierr);
6042a155e38SStefano Zampini   } else {
6052a155e38SStefano Zampini     PetscInt n_I;
6062a155e38SStefano Zampini 
60734a97f8cSStefano Zampini     /* IS for I dofs in original numbering */
60834a97f8cSStefano Zampini     ierr = PetscObjectReference((PetscObject)is_A_I);CHKERRQ(ierr);
6092a155e38SStefano Zampini     sub_schurs->is_AEj_I[0] = is_A_I;
61034a97f8cSStefano Zampini 
6112a155e38SStefano Zampini     /* IS for I dofs in I numbering (strided 1) */
6122a155e38SStefano Zampini     ierr = ISGetSize(is_A_I,&n_I);CHKERRQ(ierr);
61334a97f8cSStefano Zampini     ierr = ISCreateStride(PetscObjectComm((PetscObject)is_A_I),n_I,0,1,&is_I);CHKERRQ(ierr);
61434a97f8cSStefano Zampini 
61534a97f8cSStefano Zampini     /* II block is the same */
61634a97f8cSStefano Zampini     ierr = PetscObjectReference((PetscObject)A_II);CHKERRQ(ierr);
61734a97f8cSStefano Zampini     AE_II = A_II;
61834a97f8cSStefano Zampini   }
61934a97f8cSStefano Zampini 
6202a155e38SStefano Zampini   /* TODO: just for compatibility with the previous version, needs to be fixed */
6212a155e38SStefano Zampini   for (i=1;i<sub_schurs->n_subs;i++) {
6222a155e38SStefano Zampini     ierr = PetscObjectReference((PetscObject)sub_schurs->is_AEj_I[0]);CHKERRQ(ierr);
6232a155e38SStefano Zampini     sub_schurs->is_AEj_I[i] = sub_schurs->is_AEj_I[0];
6242a155e38SStefano Zampini   }
62534a97f8cSStefano Zampini 
6262a155e38SStefano Zampini   /* subsets in original and boundary numbering */
6272a155e38SStefano Zampini   ierr = ISLocalToGlobalMappingCreateIS(is_A_B,&BtoNmap);CHKERRQ(ierr);
6282a155e38SStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
6292a155e38SStefano Zampini     ierr = ISDuplicate(is_cc[i],&sub_schurs->is_AEj_B[i]);CHKERRQ(ierr);
6302a155e38SStefano Zampini     ierr = ISSort(sub_schurs->is_AEj_B[i]);CHKERRQ(ierr);
6312a155e38SStefano Zampini     ierr = ISGlobalToLocalMappingApplyIS(BtoNmap,IS_GTOLM_DROP,sub_schurs->is_AEj_B[i],&is_subset_B[i]);CHKERRQ(ierr);
6322a155e38SStefano Zampini   }
6332a155e38SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&BtoNmap);CHKERRQ(ierr);
6342a155e38SStefano Zampini 
6352a155e38SStefano Zampini   /* EE block */
6362a155e38SStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
6372a155e38SStefano Zampini     ierr = MatGetSubMatrix(A_BB,is_subset_B[i],is_subset_B[i],MAT_INITIAL_MATRIX,&AE_EE[i]);CHKERRQ(ierr);
6382a155e38SStefano Zampini   }
6392a155e38SStefano Zampini   /* IE block */
6402a155e38SStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
6412a155e38SStefano Zampini     ierr = MatGetSubMatrix(A_IB,is_I,is_subset_B[i],MAT_INITIAL_MATRIX,&AE_IE[i]);CHKERRQ(ierr);
6422a155e38SStefano Zampini   }
64334a97f8cSStefano Zampini   /* EI block */
6442a155e38SStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
6452a155e38SStefano Zampini     ierr = MatGetSubMatrix(A_BI,is_subset_B[i],is_I,MAT_INITIAL_MATRIX,&AE_EI[i]);CHKERRQ(ierr);
6462a155e38SStefano Zampini   }
64734a97f8cSStefano Zampini 
64834a97f8cSStefano Zampini   /* setup Schur complements on subset */
6492a155e38SStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
6502a155e38SStefano Zampini     ierr = MatCreateSchurComplement(AE_II,AE_II,AE_IE[i],AE_EI[i],AE_EE[i],&sub_schurs->S_Ej[i]);CHKERRQ(ierr);
6518a26ef87SStefano Zampini     ierr = MatCreateVecs(sub_schurs->S_Ej[i],&sub_schurs->work1[i],&sub_schurs->work2[i]);CHKERRQ(ierr);
65234a97f8cSStefano Zampini     if (AE_II == A_II) { /* we can reuse the same ksp */
65334a97f8cSStefano Zampini       KSP ksp;
65434a97f8cSStefano Zampini       ierr = MatSchurComplementGetKSP(S,&ksp);CHKERRQ(ierr);
65534a97f8cSStefano Zampini       ierr = MatSchurComplementSetKSP(sub_schurs->S_Ej[i],ksp);CHKERRQ(ierr);
65634a97f8cSStefano Zampini     } else { /* build new ksp object which inherits ksp and pc types from the original one */
65734a97f8cSStefano Zampini       KSP      origksp,schurksp;
65834a97f8cSStefano Zampini       PC       origpc,schurpc;
65934a97f8cSStefano Zampini       KSPType  ksp_type;
66034a97f8cSStefano Zampini       PCType   pc_type;
66134a97f8cSStefano Zampini       PetscInt n_internal;
66234a97f8cSStefano Zampini 
66334a97f8cSStefano Zampini       ierr = MatSchurComplementGetKSP(S,&origksp);CHKERRQ(ierr);
66434a97f8cSStefano Zampini       ierr = MatSchurComplementGetKSP(sub_schurs->S_Ej[i],&schurksp);CHKERRQ(ierr);
66534a97f8cSStefano Zampini       ierr = KSPGetType(origksp,&ksp_type);CHKERRQ(ierr);
66634a97f8cSStefano Zampini       ierr = KSPSetType(schurksp,ksp_type);CHKERRQ(ierr);
66734a97f8cSStefano Zampini       ierr = KSPGetPC(schurksp,&schurpc);CHKERRQ(ierr);
66834a97f8cSStefano Zampini       ierr = KSPGetPC(origksp,&origpc);CHKERRQ(ierr);
66934a97f8cSStefano Zampini       ierr = PCGetType(origpc,&pc_type);CHKERRQ(ierr);
67034a97f8cSStefano Zampini       ierr = PCSetType(schurpc,pc_type);CHKERRQ(ierr);
67134a97f8cSStefano Zampini       ierr = ISGetSize(is_I,&n_internal);CHKERRQ(ierr);
67234a97f8cSStefano Zampini       if (n_internal) { /* UMFPACK gives error with 0 sized problems */
67334a97f8cSStefano Zampini         MatSolverPackage solver=NULL;
67434a97f8cSStefano Zampini         ierr = PCFactorGetMatSolverPackage(origpc,(const MatSolverPackage*)&solver);CHKERRQ(ierr);
67534a97f8cSStefano Zampini         if (solver) {
67634a97f8cSStefano Zampini           ierr = PCFactorSetMatSolverPackage(schurpc,solver);CHKERRQ(ierr);
67734a97f8cSStefano Zampini         }
67834a97f8cSStefano Zampini       }
67934a97f8cSStefano Zampini       ierr = KSPSetUp(schurksp);CHKERRQ(ierr);
68034a97f8cSStefano Zampini     }
68134a97f8cSStefano Zampini   }
68234a97f8cSStefano Zampini   /* free */
6832a155e38SStefano Zampini   ierr = ISDestroy(&is_I);CHKERRQ(ierr);
6842a155e38SStefano Zampini   ierr = MatDestroy(&AE_II);CHKERRQ(ierr);
6852a155e38SStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
6862a155e38SStefano Zampini     ierr = MatDestroy(&AE_EE[i]);CHKERRQ(ierr);
6872a155e38SStefano Zampini     ierr = MatDestroy(&AE_IE[i]);CHKERRQ(ierr);
6882a155e38SStefano Zampini     ierr = MatDestroy(&AE_EI[i]);CHKERRQ(ierr);
6892a155e38SStefano Zampini     ierr = ISDestroy(&is_subset_B[i]);CHKERRQ(ierr);
6902a155e38SStefano Zampini   }
6912a155e38SStefano Zampini   ierr = PetscFree4(is_subset_B,AE_IE,AE_EI,AE_EE);CHKERRQ(ierr);
69234a97f8cSStefano Zampini   PetscFunctionReturn(0);
69334a97f8cSStefano Zampini }
694