xref: /petsc/src/ksp/pc/impls/bddc/bddcschurs.c (revision 65d8bf0a82b01a6bbcb5d0fd403f7ff82932c095)
134a97f8cSStefano Zampini #include <../src/ksp/pc/impls/bddc/bddc.h>
234a97f8cSStefano Zampini #include <../src/ksp/pc/impls/bddc/bddcprivate.h>
334a97f8cSStefano Zampini 
43202ece2SStefano Zampini PETSC_STATIC_INLINE PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt*,PetscInt,PetscBT,PetscInt*,PetscInt*,PetscInt*);
53202ece2SStefano Zampini static PetscErrorCode PCBDDCComputeExplicitSchur(Mat M, Mat *S);
63202ece2SStefano Zampini 
73202ece2SStefano Zampini #undef __FUNCT__
83202ece2SStefano Zampini #define __FUNCT__ "PCBDDCComputeExplicitSchur"
93202ece2SStefano Zampini static PetscErrorCode PCBDDCComputeExplicitSchur(Mat M, Mat *S)
103202ece2SStefano Zampini {
113202ece2SStefano Zampini   Mat            B, C, D, Bd, Cd, AinvBd;
123202ece2SStefano Zampini   KSP            ksp;
133202ece2SStefano Zampini   PC             pc;
143202ece2SStefano Zampini   PetscBool      isLU, isILU, isCHOL, Bdense, Cdense;
153202ece2SStefano Zampini   PetscReal      fill = 2.0;
163202ece2SStefano Zampini   PetscMPIInt    size;
173202ece2SStefano Zampini   PetscErrorCode ierr;
183202ece2SStefano Zampini 
193202ece2SStefano Zampini   PetscFunctionBegin;
203202ece2SStefano Zampini   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)M),&size);CHKERRQ(ierr);
213202ece2SStefano Zampini   if (size != 1) {
223202ece2SStefano Zampini     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not for parallel matrices");
233202ece2SStefano Zampini   }
243202ece2SStefano Zampini   ierr = MatSchurComplementGetSubMatrices(M, NULL, NULL, &B, &C, &D);CHKERRQ(ierr);
253202ece2SStefano Zampini   ierr = MatSchurComplementGetKSP(M, &ksp);CHKERRQ(ierr);
263202ece2SStefano Zampini   ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr);
273202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) pc, PCLU, &isLU);CHKERRQ(ierr);
283202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) pc, PCILU, &isILU);CHKERRQ(ierr);
293202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) pc, PCCHOLESKY, &isCHOL);CHKERRQ(ierr);
303202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) B, MATSEQDENSE, &Bdense);CHKERRQ(ierr);
313202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) C, MATSEQDENSE, &Cdense);CHKERRQ(ierr);
323202ece2SStefano Zampini   if (!Bdense) {
333202ece2SStefano Zampini     ierr = MatConvert(B, MATSEQDENSE, MAT_INITIAL_MATRIX, &Bd);CHKERRQ(ierr);
343202ece2SStefano Zampini   } else {
353202ece2SStefano Zampini     Bd = B;
363202ece2SStefano Zampini   }
373202ece2SStefano Zampini 
383202ece2SStefano Zampini   if (isLU || isILU || isCHOL) {
393202ece2SStefano Zampini     Mat fact;
403202ece2SStefano Zampini 
413202ece2SStefano Zampini     ierr = KSPSetUp(ksp);CHKERRQ(ierr);
423202ece2SStefano Zampini     ierr = PCFactorGetMatrix(pc, &fact);CHKERRQ(ierr);
433202ece2SStefano Zampini     ierr = MatDuplicate(Bd, MAT_DO_NOT_COPY_VALUES, &AinvBd);CHKERRQ(ierr);
443202ece2SStefano Zampini     ierr = MatMatSolve(fact, Bd, AinvBd);CHKERRQ(ierr);
453202ece2SStefano Zampini   } else {
463202ece2SStefano Zampini     Mat Ainvd;
473202ece2SStefano Zampini 
483202ece2SStefano Zampini     ierr = PCComputeExplicitOperator(pc, &Ainvd);CHKERRQ(ierr);
493202ece2SStefano Zampini     ierr = MatMatMult(Ainvd, Bd, MAT_INITIAL_MATRIX, fill, &AinvBd);CHKERRQ(ierr);
503202ece2SStefano Zampini     ierr = MatDestroy(&Ainvd);CHKERRQ(ierr);
513202ece2SStefano Zampini   }
523202ece2SStefano Zampini   if (!Bdense) {
533202ece2SStefano Zampini     ierr = MatDestroy(&Bd);CHKERRQ(ierr);
543202ece2SStefano Zampini   }
553202ece2SStefano Zampini   if (!Cdense) {
563202ece2SStefano Zampini     ierr = MatConvert(C, MATSEQDENSE, MAT_INITIAL_MATRIX, &Cd);CHKERRQ(ierr);
573202ece2SStefano Zampini   } else {
583202ece2SStefano Zampini     Cd = C;
593202ece2SStefano Zampini   }
603202ece2SStefano Zampini 
613202ece2SStefano Zampini   ierr = MatMatMult(Cd, AinvBd, MAT_INITIAL_MATRIX, fill, S);CHKERRQ(ierr);
623202ece2SStefano Zampini   ierr = MatDestroy(&AinvBd);CHKERRQ(ierr);
633202ece2SStefano Zampini   if (!Cdense) {
643202ece2SStefano Zampini     ierr = MatDestroy(&Cd);CHKERRQ(ierr);
653202ece2SStefano Zampini   }
663202ece2SStefano Zampini 
673202ece2SStefano Zampini   if (D) {
683202ece2SStefano Zampini     Mat       Dd;
693202ece2SStefano Zampini     PetscBool Ddense;
703202ece2SStefano Zampini 
713202ece2SStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)D,MATSEQDENSE,&Ddense);CHKERRQ(ierr);
723202ece2SStefano Zampini     if (!Ddense) {
733202ece2SStefano Zampini       ierr = MatConvert(D, MATSEQDENSE, MAT_INITIAL_MATRIX, &Dd);CHKERRQ(ierr);
743202ece2SStefano Zampini     } else {
753202ece2SStefano Zampini       Dd = D;
763202ece2SStefano Zampini     }
773202ece2SStefano Zampini     ierr = MatAYPX(*S,-1.0,Dd,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
783202ece2SStefano Zampini     if (!Ddense) {
793202ece2SStefano Zampini       ierr = MatDestroy(&Dd);CHKERRQ(ierr);
803202ece2SStefano Zampini     }
813202ece2SStefano Zampini   } else {
823202ece2SStefano Zampini     ierr = MatScale(*S,-1.0);CHKERRQ(ierr);
833202ece2SStefano Zampini   }
843202ece2SStefano Zampini   PetscFunctionReturn(0);
853202ece2SStefano Zampini }
8634a97f8cSStefano Zampini 
8734a97f8cSStefano Zampini #undef __FUNCT__
881580ed26SStefano Zampini #define __FUNCT__ "PCBDDCSubSchursSetUp"
891580ed26SStefano Zampini PetscErrorCode PCBDDCSubSchursSetUp(PCBDDCSubSchurs sub_schurs, PetscInt xadj[], PetscInt adjncy[], PetscInt nlayers)
90b1b3d7a2SStefano Zampini {
91b1b3d7a2SStefano Zampini   Mat                    A_II,A_IB,A_BI,A_BB;
92b1b3d7a2SStefano Zampini   Mat                    AE_II,*AE_IE,*AE_EI,*AE_EE;
93883469d8SStefano Zampini   Mat                    S_all,global_schur_subsets,work_mat;
945db18549SStefano Zampini   ISLocalToGlobalMapping l2gmap_subsets;
955db18549SStefano Zampini   IS                     is_I,*is_subset_B,temp_is;
96883469d8SStefano Zampini   PetscInt               *nnz,*all_local_idx_G,*all_local_idx_B,*all_local_idx_N;
975db18549SStefano Zampini   PetscInt               i,subset_size,max_subset_size;
98883469d8SStefano Zampini   PetscInt               extra,local_size,global_size;
99b1b3d7a2SStefano Zampini   PetscErrorCode         ierr;
100b1b3d7a2SStefano Zampini 
101b1b3d7a2SStefano Zampini   PetscFunctionBegin;
102883469d8SStefano Zampini 
103b1b3d7a2SStefano Zampini   /* allocate space for schur complements */
10468270318SStefano Zampini   ierr = PetscMalloc4(sub_schurs->n_subs,&sub_schurs->is_AEj_B,
105b1b3d7a2SStefano Zampini                       sub_schurs->n_subs,&sub_schurs->S_Ej,
106b1b3d7a2SStefano Zampini                       sub_schurs->n_subs,&sub_schurs->work1,
107b1b3d7a2SStefano Zampini                       sub_schurs->n_subs,&sub_schurs->work2);CHKERRQ(ierr);
108b1b3d7a2SStefano Zampini 
109b1b3d7a2SStefano Zampini   /* get Schur complement matrices */
110883469d8SStefano Zampini   if (!sub_schurs->use_mumps) {
111b1b3d7a2SStefano Zampini     ierr = MatSchurComplementGetSubMatrices(sub_schurs->S,&A_II,NULL,&A_IB,&A_BI,&A_BB);CHKERRQ(ierr);
112b1b3d7a2SStefano Zampini     ierr = PetscMalloc4(sub_schurs->n_subs,&is_subset_B,
113b1b3d7a2SStefano Zampini                         sub_schurs->n_subs,&AE_IE,
114b1b3d7a2SStefano Zampini                         sub_schurs->n_subs,&AE_EI,
115b1b3d7a2SStefano Zampini                         sub_schurs->n_subs,&AE_EE);CHKERRQ(ierr);
116b1b3d7a2SStefano Zampini   }
117b1b3d7a2SStefano Zampini 
118b1b3d7a2SStefano Zampini   /* determine interior problems */
119b1b3d7a2SStefano Zampini   if (nlayers >= 0 && xadj != NULL && adjncy != NULL) { /* Interior problems can be different from the original one */
120b1b3d7a2SStefano Zampini     PetscBT                touched;
121b1b3d7a2SStefano Zampini     const PetscInt*        idx_B;
122b1b3d7a2SStefano Zampini     PetscInt               n_I,n_B,n_local_dofs,n_prev_added,j,layer,*local_numbering;
123b1b3d7a2SStefano Zampini 
124b1b3d7a2SStefano Zampini     /* get sizes */
125b1b3d7a2SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr);
126b1b3d7a2SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_B,&n_B);CHKERRQ(ierr);
127b1b3d7a2SStefano Zampini 
128b1b3d7a2SStefano Zampini     ierr = PetscMalloc1(n_I+n_B,&local_numbering);CHKERRQ(ierr);
129b1b3d7a2SStefano Zampini     ierr = PetscBTCreate(n_I+n_B,&touched);CHKERRQ(ierr);
130b1b3d7a2SStefano Zampini     ierr = PetscBTMemzero(n_I+n_B,touched);CHKERRQ(ierr);
131b1b3d7a2SStefano Zampini 
132b1b3d7a2SStefano Zampini     /* all boundary dofs must be skipped when adding layers */
133b1b3d7a2SStefano Zampini     ierr = ISGetIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr);
134b1b3d7a2SStefano Zampini     for (j=0;j<n_B;j++) {
135b1b3d7a2SStefano Zampini       ierr = PetscBTSet(touched,idx_B[j]);CHKERRQ(ierr);
136b1b3d7a2SStefano Zampini     }
137b1b3d7a2SStefano Zampini     ierr = PetscMemcpy(local_numbering,idx_B,n_B*sizeof(PetscInt));CHKERRQ(ierr);
138b1b3d7a2SStefano Zampini     ierr = ISRestoreIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr);
139b1b3d7a2SStefano Zampini 
140b1b3d7a2SStefano Zampini     /* add prescribed number of layers of dofs */
141b1b3d7a2SStefano Zampini     n_local_dofs = n_B;
142b1b3d7a2SStefano Zampini     n_prev_added = n_B;
143b1b3d7a2SStefano Zampini     for (layer=0;layer<nlayers;layer++) {
144b1b3d7a2SStefano Zampini       PetscInt n_added;
145b1b3d7a2SStefano Zampini       if (n_local_dofs == n_I+n_B) break;
146b1b3d7a2SStefano Zampini       if (n_local_dofs > n_I+n_B) {
147b1b3d7a2SStefano 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);
148b1b3d7a2SStefano Zampini       }
149b1b3d7a2SStefano Zampini       ierr = PCBDDCAdjGetNextLayer_Private(local_numbering+n_local_dofs,n_prev_added,touched,xadj,adjncy,&n_added);CHKERRQ(ierr);
150b1b3d7a2SStefano Zampini       n_prev_added = n_added;
151b1b3d7a2SStefano Zampini       n_local_dofs += n_added;
152b1b3d7a2SStefano Zampini       if (!n_added) break;
153b1b3d7a2SStefano Zampini     }
154b1b3d7a2SStefano Zampini     ierr = PetscBTDestroy(&touched);CHKERRQ(ierr);
155b1b3d7a2SStefano Zampini 
156883469d8SStefano Zampini     /* IS for I layer dofs in original numbering */
15768270318SStefano 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);
158b1b3d7a2SStefano Zampini     ierr = PetscFree(local_numbering);CHKERRQ(ierr);
15968270318SStefano Zampini     ierr = ISSort(sub_schurs->is_I_layer);CHKERRQ(ierr);
160883469d8SStefano Zampini     /* IS for I layer dofs in I numbering */
161883469d8SStefano Zampini     if (!sub_schurs->use_mumps) {
162b1b3d7a2SStefano Zampini       ISLocalToGlobalMapping ItoNmap;
163b1b3d7a2SStefano Zampini       ierr = ISLocalToGlobalMappingCreateIS(sub_schurs->is_I,&ItoNmap);CHKERRQ(ierr);
16468270318SStefano Zampini       ierr = ISGlobalToLocalMappingApplyIS(ItoNmap,IS_GTOLM_DROP,sub_schurs->is_I_layer,&is_I);CHKERRQ(ierr);
165b1b3d7a2SStefano Zampini       ierr = ISLocalToGlobalMappingDestroy(&ItoNmap);CHKERRQ(ierr);
166b1b3d7a2SStefano Zampini 
167b1b3d7a2SStefano Zampini       /* II block */
168b1b3d7a2SStefano Zampini       ierr = MatGetSubMatrix(A_II,is_I,is_I,MAT_INITIAL_MATRIX,&AE_II);CHKERRQ(ierr);
169b1b3d7a2SStefano Zampini     }
170b1b3d7a2SStefano Zampini   } else {
171b1b3d7a2SStefano Zampini     PetscInt n_I;
172b1b3d7a2SStefano Zampini 
173b1b3d7a2SStefano Zampini     /* IS for I dofs in original numbering */
174b1b3d7a2SStefano Zampini     ierr = PetscObjectReference((PetscObject)sub_schurs->is_I);CHKERRQ(ierr);
17568270318SStefano Zampini     sub_schurs->is_I_layer = sub_schurs->is_I;
176b1b3d7a2SStefano Zampini 
177b1b3d7a2SStefano Zampini     /* IS for I dofs in I numbering (strided 1) */
178883469d8SStefano Zampini     if (!sub_schurs->use_mumps) {
179b1b3d7a2SStefano Zampini       ierr = ISGetSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr);
180b1b3d7a2SStefano Zampini       ierr = ISCreateStride(PetscObjectComm((PetscObject)sub_schurs->is_I),n_I,0,1,&is_I);CHKERRQ(ierr);
181b1b3d7a2SStefano Zampini 
182b1b3d7a2SStefano Zampini       /* II block is the same */
183b1b3d7a2SStefano Zampini       ierr = PetscObjectReference((PetscObject)A_II);CHKERRQ(ierr);
184b1b3d7a2SStefano Zampini       AE_II = A_II;
185b1b3d7a2SStefano Zampini     }
186b1b3d7a2SStefano Zampini   }
187b1b3d7a2SStefano Zampini 
188b1b3d7a2SStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
189b1b3d7a2SStefano Zampini     ierr = ISDuplicate(sub_schurs->is_subs[i],&sub_schurs->is_AEj_B[i]);CHKERRQ(ierr);
190b1b3d7a2SStefano Zampini     ierr = ISSort(sub_schurs->is_AEj_B[i]);CHKERRQ(ierr);
191883469d8SStefano Zampini   }
192883469d8SStefano Zampini 
193883469d8SStefano Zampini   /* Get info on subset sizes and sum of all subsets sizes */
194883469d8SStefano Zampini   max_subset_size = 0;
195883469d8SStefano Zampini   local_size = 0;
196883469d8SStefano Zampini   for (i=0;i<sub_schurs->n_subs_seq;i++) {
197883469d8SStefano Zampini     PetscInt j = sub_schurs->index_sequential[i];
198883469d8SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_AEj_B[j],&subset_size);CHKERRQ(ierr);
199883469d8SStefano Zampini     max_subset_size = PetscMax(subset_size,max_subset_size);
200883469d8SStefano Zampini     local_size += subset_size;
201883469d8SStefano Zampini   }
202883469d8SStefano Zampini 
203883469d8SStefano Zampini   /* Work arrays for local indices */
204883469d8SStefano Zampini   ierr = PetscMalloc1(local_size,&all_local_idx_B);CHKERRQ(ierr);
205883469d8SStefano Zampini   extra = 0;
206883469d8SStefano Zampini   if (sub_schurs->use_mumps) {
207883469d8SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_I_layer,&extra);CHKERRQ(ierr);
208883469d8SStefano Zampini   }
209883469d8SStefano Zampini   ierr = PetscMalloc1(local_size+extra,&all_local_idx_N);CHKERRQ(ierr);
210883469d8SStefano Zampini   if (extra) {
211883469d8SStefano Zampini     const PetscInt *idxs;
212883469d8SStefano Zampini     ierr = ISGetIndices(sub_schurs->is_I_layer,&idxs);CHKERRQ(ierr);
213883469d8SStefano Zampini     ierr = PetscMemcpy(all_local_idx_N,idxs,extra*sizeof(PetscInt));CHKERRQ(ierr);
214883469d8SStefano Zampini     ierr = ISRestoreIndices(sub_schurs->is_I_layer,&idxs);CHKERRQ(ierr);
215883469d8SStefano Zampini   }
216883469d8SStefano Zampini   ierr = PetscMalloc1(local_size,&nnz);CHKERRQ(ierr);
217883469d8SStefano Zampini 
218883469d8SStefano Zampini   /* Get local indices in local numbering */
219883469d8SStefano Zampini   local_size = 0;
220883469d8SStefano Zampini   for (i=0;i<sub_schurs->n_subs_seq;i++) {
221883469d8SStefano Zampini     PetscInt j;
222883469d8SStefano Zampini     const    PetscInt *idxs;
223883469d8SStefano Zampini 
224883469d8SStefano Zampini     PetscInt local_problem_index = sub_schurs->index_sequential[i];
225883469d8SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_AEj_B[local_problem_index],&subset_size);CHKERRQ(ierr);
226883469d8SStefano Zampini     ierr = ISGetIndices(sub_schurs->is_AEj_B[local_problem_index],&idxs);CHKERRQ(ierr);
227883469d8SStefano Zampini     /* subset indices in local numbering */
228883469d8SStefano Zampini     ierr = PetscMemcpy(all_local_idx_N+local_size+extra,idxs,subset_size*sizeof(PetscInt));CHKERRQ(ierr);
229883469d8SStefano Zampini     ierr = ISRestoreIndices(sub_schurs->is_AEj_B[local_problem_index],&idxs);CHKERRQ(ierr);
230883469d8SStefano Zampini     for (j=0;j<subset_size;j++) nnz[local_size+j] = subset_size;
231883469d8SStefano Zampini     local_size += subset_size;
232883469d8SStefano Zampini   }
233883469d8SStefano Zampini 
234883469d8SStefano Zampini   S_all = NULL;
235883469d8SStefano Zampini   if (!sub_schurs->use_mumps) {
236883469d8SStefano Zampini     /* subsets in original and boundary numbering */
237883469d8SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
2385db18549SStefano Zampini       ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,sub_schurs->is_AEj_B[i],&is_subset_B[i]);CHKERRQ(ierr);
239b1b3d7a2SStefano Zampini     }
240b1b3d7a2SStefano Zampini 
241b1b3d7a2SStefano Zampini     /* EE block */
242b1b3d7a2SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
243b1b3d7a2SStefano Zampini       ierr = MatGetSubMatrix(A_BB,is_subset_B[i],is_subset_B[i],MAT_INITIAL_MATRIX,&AE_EE[i]);CHKERRQ(ierr);
244b1b3d7a2SStefano Zampini     }
245b1b3d7a2SStefano Zampini     /* IE block */
246b1b3d7a2SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
247b1b3d7a2SStefano Zampini       ierr = MatGetSubMatrix(A_IB,is_I,is_subset_B[i],MAT_INITIAL_MATRIX,&AE_IE[i]);CHKERRQ(ierr);
248b1b3d7a2SStefano Zampini     }
249b1b3d7a2SStefano Zampini     /* EI block */
250b1b3d7a2SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
251b1b3d7a2SStefano Zampini       ierr = MatGetSubMatrix(A_BI,is_subset_B[i],is_I,MAT_INITIAL_MATRIX,&AE_EI[i]);CHKERRQ(ierr);
252b1b3d7a2SStefano Zampini     }
253b1b3d7a2SStefano Zampini 
254b1b3d7a2SStefano Zampini     /* setup Schur complements on subset */
255b1b3d7a2SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
256b1b3d7a2SStefano Zampini       ierr = MatCreateSchurComplement(AE_II,AE_II,AE_IE[i],AE_EI[i],AE_EE[i],&sub_schurs->S_Ej[i]);CHKERRQ(ierr);
257b1b3d7a2SStefano Zampini       if (AE_II == A_II) { /* we can reuse the same ksp */
258b1b3d7a2SStefano Zampini         KSP ksp;
259b1b3d7a2SStefano Zampini         ierr = MatSchurComplementGetKSP(sub_schurs->S,&ksp);CHKERRQ(ierr);
260b1b3d7a2SStefano Zampini         ierr = MatSchurComplementSetKSP(sub_schurs->S_Ej[i],ksp);CHKERRQ(ierr);
261b1b3d7a2SStefano Zampini       } else { /* build new ksp object which inherits ksp and pc types from the original one */
262b1b3d7a2SStefano Zampini         KSP      origksp,schurksp;
263b1b3d7a2SStefano Zampini         PC       origpc,schurpc;
264b1b3d7a2SStefano Zampini         KSPType  ksp_type;
265b1b3d7a2SStefano Zampini         PCType   pc_type;
266b1b3d7a2SStefano Zampini         PetscInt n_internal;
267b1b3d7a2SStefano Zampini 
268b1b3d7a2SStefano Zampini         ierr = MatSchurComplementGetKSP(sub_schurs->S,&origksp);CHKERRQ(ierr);
269b1b3d7a2SStefano Zampini         ierr = MatSchurComplementGetKSP(sub_schurs->S_Ej[i],&schurksp);CHKERRQ(ierr);
270b1b3d7a2SStefano Zampini         ierr = KSPGetType(origksp,&ksp_type);CHKERRQ(ierr);
271b1b3d7a2SStefano Zampini         ierr = KSPSetType(schurksp,ksp_type);CHKERRQ(ierr);
272b1b3d7a2SStefano Zampini         ierr = KSPGetPC(schurksp,&schurpc);CHKERRQ(ierr);
273b1b3d7a2SStefano Zampini         ierr = KSPGetPC(origksp,&origpc);CHKERRQ(ierr);
274b1b3d7a2SStefano Zampini         ierr = PCGetType(origpc,&pc_type);CHKERRQ(ierr);
275b1b3d7a2SStefano Zampini         ierr = PCSetType(schurpc,pc_type);CHKERRQ(ierr);
276b1b3d7a2SStefano Zampini         ierr = ISGetSize(is_I,&n_internal);CHKERRQ(ierr);
277b1b3d7a2SStefano Zampini         if (n_internal) { /* UMFPACK gives error with 0 sized problems */
278b1b3d7a2SStefano Zampini           MatSolverPackage solver=NULL;
279b1b3d7a2SStefano Zampini           ierr = PCFactorGetMatSolverPackage(origpc,(const MatSolverPackage*)&solver);CHKERRQ(ierr);
280b1b3d7a2SStefano Zampini           if (solver) {
281b1b3d7a2SStefano Zampini             ierr = PCFactorSetMatSolverPackage(schurpc,solver);CHKERRQ(ierr);
282b1b3d7a2SStefano Zampini           }
283b1b3d7a2SStefano Zampini         }
284b1b3d7a2SStefano Zampini         ierr = KSPSetUp(schurksp);CHKERRQ(ierr);
285b1b3d7a2SStefano Zampini       }
286b1b3d7a2SStefano Zampini     }
287b1b3d7a2SStefano Zampini     /* free */
288b1b3d7a2SStefano Zampini     ierr = ISDestroy(&is_I);CHKERRQ(ierr);
289b1b3d7a2SStefano Zampini     ierr = MatDestroy(&AE_II);CHKERRQ(ierr);
290b1b3d7a2SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
291b1b3d7a2SStefano Zampini       ierr = MatDestroy(&AE_EE[i]);CHKERRQ(ierr);
292b1b3d7a2SStefano Zampini       ierr = MatDestroy(&AE_IE[i]);CHKERRQ(ierr);
293b1b3d7a2SStefano Zampini       ierr = MatDestroy(&AE_EI[i]);CHKERRQ(ierr);
294b1b3d7a2SStefano Zampini       ierr = ISDestroy(&is_subset_B[i]);CHKERRQ(ierr);
295b1b3d7a2SStefano Zampini     }
296b1b3d7a2SStefano Zampini     ierr = PetscFree4(is_subset_B,AE_IE,AE_EI,AE_EE);CHKERRQ(ierr);
297883469d8SStefano Zampini #if defined(PETSC_HAVE_MUMPS)
298883469d8SStefano Zampini   } else {
299883469d8SStefano Zampini     Mat           A,F;
300883469d8SStefano Zampini     IS            is_A_all;
301883469d8SStefano Zampini     PetscBool     is_symmetric;
302883469d8SStefano Zampini     PetscInt      *idxs_schur,n_I;
303883469d8SStefano Zampini 
304883469d8SStefano Zampini     /* get working mat */
305883469d8SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_I_layer,&n_I);CHKERRQ(ierr);
306883469d8SStefano Zampini     ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size+n_I,all_local_idx_N,PETSC_COPY_VALUES,&is_A_all);CHKERRQ(ierr);
307883469d8SStefano Zampini     ierr = MatGetSubMatrixUnsorted(sub_schurs->A,is_A_all,is_A_all,&A);CHKERRQ(ierr);
308883469d8SStefano Zampini     ierr = ISDestroy(&is_A_all);CHKERRQ(ierr);
309883469d8SStefano Zampini 
310883469d8SStefano Zampini     ierr = MatIsSymmetric(sub_schurs->A,0.0,&is_symmetric);CHKERRQ(ierr);
311883469d8SStefano Zampini     if (is_symmetric) {
312883469d8SStefano Zampini       ierr = MatGetFactor(A,MATSOLVERMUMPS,MAT_FACTOR_CHOLESKY,&F);CHKERRQ(ierr);
313883469d8SStefano Zampini     } else {
314883469d8SStefano Zampini       ierr = MatGetFactor(A,MATSOLVERMUMPS,MAT_FACTOR_LU,&F);CHKERRQ(ierr);
315883469d8SStefano Zampini     }
316883469d8SStefano Zampini 
317883469d8SStefano Zampini     /* subsets ordered last */
318883469d8SStefano Zampini     ierr = PetscMalloc1(local_size,&idxs_schur);CHKERRQ(ierr);
319883469d8SStefano Zampini     for (i=0;i<local_size;i++) {
320883469d8SStefano Zampini       idxs_schur[i] = n_I+i+1;
321883469d8SStefano Zampini     }
322883469d8SStefano Zampini     ierr = MatMumpsSetSchurIndices(F,local_size,idxs_schur);CHKERRQ(ierr);
323883469d8SStefano Zampini     ierr = PetscFree(idxs_schur);CHKERRQ(ierr);
324883469d8SStefano Zampini 
325883469d8SStefano Zampini     /* factorization step */
326883469d8SStefano Zampini     if (is_symmetric) {
327883469d8SStefano Zampini       ierr = MatCholeskyFactorSymbolic(F,A,NULL,NULL);CHKERRQ(ierr);
328883469d8SStefano Zampini       ierr = MatCholeskyFactorNumeric(F,A,NULL);CHKERRQ(ierr);
329883469d8SStefano Zampini     } else {
330883469d8SStefano Zampini       ierr = MatLUFactorSymbolic(F,A,NULL,NULL,NULL);CHKERRQ(ierr);
331883469d8SStefano Zampini       ierr = MatLUFactorNumeric(F,A,NULL);CHKERRQ(ierr);
332883469d8SStefano Zampini     }
333883469d8SStefano Zampini 
334883469d8SStefano Zampini     /* get explicit Schur Complement computed during numeric factorization */
335883469d8SStefano Zampini     ierr = MatMumpsGetSchurComplement(F,&S_all);CHKERRQ(ierr);
336883469d8SStefano Zampini 
337883469d8SStefano Zampini     /* free workspace */
338883469d8SStefano Zampini     ierr = MatDestroy(&A);CHKERRQ(ierr);
339883469d8SStefano Zampini     ierr = MatDestroy(&F);CHKERRQ(ierr);
340883469d8SStefano Zampini 
341883469d8SStefano Zampini #endif
342b1b3d7a2SStefano Zampini   }
3435db18549SStefano Zampini 
3449221af80SStefano Zampini   /* TODO: just for compatibility with the previous version, needs to be fixed */
3459221af80SStefano Zampini   for (i=0;i<sub_schurs->n_subs_par;i++) {
3469221af80SStefano Zampini     PetscInt j = sub_schurs->index_parallel[i];
3479221af80SStefano Zampini     ierr = MatCreateVecs(sub_schurs->S_Ej[j],&sub_schurs->work1[j],&sub_schurs->work2[j]);CHKERRQ(ierr);
3489221af80SStefano Zampini   }
3499221af80SStefano Zampini   for (i=0;i<sub_schurs->n_subs_seq;i++) {
3509221af80SStefano Zampini      sub_schurs->work1[sub_schurs->index_sequential[i]] = 0;
3519221af80SStefano Zampini      sub_schurs->work2[sub_schurs->index_sequential[i]] = 0;
3529221af80SStefano Zampini   }
3539221af80SStefano Zampini 
3545db18549SStefano Zampini   if (!sub_schurs->n_subs_seq_g) {
3559221af80SStefano Zampini     sub_schurs->S_Ej_all = 0;
3569221af80SStefano Zampini     sub_schurs->sum_S_Ej_all = 0;
3579221af80SStefano Zampini     sub_schurs->is_Ej_all = 0;
3585db18549SStefano Zampini     PetscFunctionReturn(0);
3595db18549SStefano Zampini   }
3605db18549SStefano Zampini 
3615db18549SStefano Zampini   /* subset indices in local boundary numbering */
362883469d8SStefano Zampini   ierr = ISGlobalToLocalMappingApply(sub_schurs->BtoNmap,IS_GTOLM_DROP,local_size,all_local_idx_N+extra,&subset_size,all_local_idx_B);CHKERRQ(ierr);
3635db18549SStefano Zampini   if (subset_size != local_size) {
3645db18549SStefano Zampini     SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in sub_schurs serial (BtoNmap)! %d != %d\n",subset_size,local_size);
3655db18549SStefano Zampini   }
3665db18549SStefano Zampini   ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size,all_local_idx_B,PETSC_OWN_POINTER,&sub_schurs->is_Ej_all);CHKERRQ(ierr);
3675db18549SStefano Zampini 
3683202ece2SStefano Zampini   /* Local matrix of all local Schur on subsets transposed */
3695db18549SStefano Zampini   ierr = MatCreate(PETSC_COMM_SELF,&sub_schurs->S_Ej_all);CHKERRQ(ierr);
3705db18549SStefano Zampini   ierr = MatSetSizes(sub_schurs->S_Ej_all,PETSC_DECIDE,PETSC_DECIDE,local_size,local_size);CHKERRQ(ierr);
3715db18549SStefano Zampini   ierr = MatSetType(sub_schurs->S_Ej_all,MATAIJ);CHKERRQ(ierr);
3725db18549SStefano Zampini   ierr = MatSeqAIJSetPreallocation(sub_schurs->S_Ej_all,0,nnz);CHKERRQ(ierr);
3735db18549SStefano Zampini   ierr = PetscFree(nnz);CHKERRQ(ierr);
3745db18549SStefano Zampini 
375883469d8SStefano Zampini   if (!sub_schurs->use_mumps) {
3765db18549SStefano Zampini     PetscScalar *fill_vals;
3775db18549SStefano Zampini     PetscInt    *dummy_idx;
3785db18549SStefano Zampini 
3795db18549SStefano Zampini     /* Work arrays */
3803202ece2SStefano Zampini     ierr = PetscMalloc1(max_subset_size,&dummy_idx);CHKERRQ(ierr);
3815db18549SStefano Zampini 
3823202ece2SStefano Zampini     /* Loop on local problems end compute Schur complements explicitly */
3835db18549SStefano Zampini     local_size = 0;
3845db18549SStefano Zampini     for (i=0;i<sub_schurs->n_subs_seq;i++) {
3853202ece2SStefano Zampini       Mat       S_Ej_expl;
3863202ece2SStefano Zampini       PetscInt  j,lpi = sub_schurs->index_sequential[i];
3873202ece2SStefano Zampini       PetscBool Sdense;
3885db18549SStefano Zampini 
3893202ece2SStefano Zampini       ierr = PCBDDCComputeExplicitSchur(sub_schurs->S_Ej[lpi],&S_Ej_expl);CHKERRQ(ierr);
3903202ece2SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_AEj_B[lpi],&subset_size);CHKERRQ(ierr);
3913202ece2SStefano Zampini       ierr = PetscObjectTypeCompare((PetscObject)S_Ej_expl,MATSEQDENSE,&Sdense);CHKERRQ(ierr);
3923202ece2SStefano Zampini       if (Sdense) {
3935db18549SStefano Zampini         for (j=0;j<subset_size;j++) {
3945db18549SStefano Zampini           dummy_idx[j]=local_size+j;
3955db18549SStefano Zampini         }
3963202ece2SStefano Zampini         ierr = MatDenseGetArray(S_Ej_expl,&fill_vals);CHKERRQ(ierr);
3975db18549SStefano Zampini         ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,fill_vals,INSERT_VALUES);CHKERRQ(ierr);
3983202ece2SStefano Zampini         ierr = MatDenseRestoreArray(S_Ej_expl,&fill_vals);CHKERRQ(ierr);
3993202ece2SStefano Zampini       } else {
4003202ece2SStefano Zampini         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not yet implemented for sparse matrices");
4015db18549SStefano Zampini       }
4023202ece2SStefano Zampini       local_size += subset_size;
403*65d8bf0aSStefano Zampini       ierr = MatDestroy(&sub_schurs->S_Ej[lpi]);CHKERRQ(ierr);
404*65d8bf0aSStefano Zampini       sub_schurs->S_Ej[lpi] = S_Ej_expl;
4053202ece2SStefano Zampini     }
4063202ece2SStefano Zampini     ierr = PetscFree(dummy_idx);CHKERRQ(ierr);
407883469d8SStefano Zampini   } else {
408883469d8SStefano Zampini     PetscInt  *dummy_idx,n_all;
409*65d8bf0aSStefano Zampini     PetscBool is_symmetric,compute_Stilda=PETSC_TRUE;
410883469d8SStefano Zampini 
411883469d8SStefano Zampini     /* Work arrays */
412883469d8SStefano Zampini     ierr = PetscMalloc1(max_subset_size,&dummy_idx);CHKERRQ(ierr);
413883469d8SStefano Zampini 
414*65d8bf0aSStefano Zampini     /* Work arrays */
415*65d8bf0aSStefano Zampini     ierr = MatGetSize(S_all,&n_all,NULL);CHKERRQ(ierr);
416*65d8bf0aSStefano Zampini     local_size = 0;
417*65d8bf0aSStefano Zampini     ierr = MatIsSymmetric(sub_schurs->A,0.0,&is_symmetric);CHKERRQ(ierr);
418*65d8bf0aSStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
419*65d8bf0aSStefano Zampini       Mat S_EE;
420*65d8bf0aSStefano Zampini       IS  is_E;
421*65d8bf0aSStefano Zampini       PetscScalar *vals;
422*65d8bf0aSStefano Zampini       PetscInt j;
423*65d8bf0aSStefano Zampini 
424*65d8bf0aSStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_AEj_B[i],&subset_size);CHKERRQ(ierr);
425*65d8bf0aSStefano Zampini       ierr = ISCreateStride(PetscObjectComm((PetscObject)sub_schurs->is_I),subset_size,local_size,1,&is_E);CHKERRQ(ierr);
426*65d8bf0aSStefano Zampini       ierr = MatGetSubMatrix(S_all,is_E,is_E,MAT_INITIAL_MATRIX,&S_EE);CHKERRQ(ierr);
427*65d8bf0aSStefano Zampini       sub_schurs->S_Ej[i] = S_EE;
428*65d8bf0aSStefano Zampini       if (compute_Stilda) {
429*65d8bf0aSStefano Zampini         Mat Stilda;
430*65d8bf0aSStefano Zampini         if (sub_schurs->n_subs == 1) {
431*65d8bf0aSStefano Zampini           ierr = PetscObjectReference((PetscObject)sub_schurs->S_Ej[i]);CHKERRQ(ierr);
432*65d8bf0aSStefano Zampini           Stilda = sub_schurs->S_Ej[i];
433*65d8bf0aSStefano Zampini         } else {
434*65d8bf0aSStefano Zampini           KSP ksp;
435*65d8bf0aSStefano Zampini           PC  pc;
436*65d8bf0aSStefano Zampini           Mat S_EF,S_FE,S_FF,Stilda_impl;
437*65d8bf0aSStefano Zampini           IS  is_F;
438*65d8bf0aSStefano Zampini 
439*65d8bf0aSStefano Zampini           ierr = ISComplement(is_E,0,n_all,&is_F);CHKERRQ(ierr);
440*65d8bf0aSStefano Zampini           ierr = MatGetSubMatrix(S_all,is_E,is_F,MAT_INITIAL_MATRIX,&S_EF);CHKERRQ(ierr);
441*65d8bf0aSStefano Zampini           ierr = MatGetSubMatrix(S_all,is_F,is_F,MAT_INITIAL_MATRIX,&S_FF);CHKERRQ(ierr);
442*65d8bf0aSStefano Zampini           ierr = MatGetSubMatrix(S_all,is_F,is_E,MAT_INITIAL_MATRIX,&S_FE);CHKERRQ(ierr);
443*65d8bf0aSStefano Zampini           ierr = MatCreateSchurComplement(S_FF,S_FF,S_FE,S_EF,S_EE,&Stilda_impl);CHKERRQ(ierr);
444*65d8bf0aSStefano Zampini           ierr = MatSchurComplementGetKSP(Stilda_impl,&ksp);CHKERRQ(ierr);
445*65d8bf0aSStefano Zampini           ierr = KSPSetType(ksp,KSPPREONLY);CHKERRQ(ierr);
446*65d8bf0aSStefano Zampini           ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr);
447*65d8bf0aSStefano Zampini           if (is_symmetric) {
448*65d8bf0aSStefano Zampini             ierr = PCSetType(pc,PCCHOLESKY);CHKERRQ(ierr);
449*65d8bf0aSStefano Zampini           } else {
450*65d8bf0aSStefano Zampini             ierr = PCSetType(pc,PCLU);CHKERRQ(ierr);
451*65d8bf0aSStefano Zampini           }
452*65d8bf0aSStefano Zampini           ierr = PCFactorSetUseInPlace(pc,PETSC_TRUE);CHKERRQ(ierr);
453*65d8bf0aSStefano Zampini           ierr = KSPSetUp(ksp);CHKERRQ(ierr);
454*65d8bf0aSStefano Zampini           ierr = PCBDDCComputeExplicitSchur(Stilda_impl,&Stilda);CHKERRQ(ierr);
455*65d8bf0aSStefano Zampini           ierr = MatDestroy(&S_FF);CHKERRQ(ierr);
456*65d8bf0aSStefano Zampini           ierr = MatDestroy(&S_FE);CHKERRQ(ierr);
457*65d8bf0aSStefano Zampini           ierr = MatDestroy(&S_EF);CHKERRQ(ierr);
458*65d8bf0aSStefano Zampini           ierr = MatDestroy(&Stilda_impl);CHKERRQ(ierr);
459*65d8bf0aSStefano Zampini           ierr = ISDestroy(&is_F);CHKERRQ(ierr);
460*65d8bf0aSStefano Zampini         }
461*65d8bf0aSStefano Zampini         ierr = MatDestroy(&Stilda);CHKERRQ(ierr);
462*65d8bf0aSStefano Zampini       }
463*65d8bf0aSStefano Zampini 
464883469d8SStefano Zampini       for (j=0;j<subset_size;j++) {
465883469d8SStefano Zampini         dummy_idx[j]=local_size+j;
466883469d8SStefano Zampini       }
467*65d8bf0aSStefano Zampini       ierr = MatDenseGetArray(S_EE,&vals);CHKERRQ(ierr);
468*65d8bf0aSStefano Zampini       ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,vals,INSERT_VALUES);CHKERRQ(ierr);
469*65d8bf0aSStefano Zampini       ierr = MatDenseRestoreArray(S_EE,&vals);CHKERRQ(ierr);
470*65d8bf0aSStefano Zampini       ierr = ISDestroy(&is_E);CHKERRQ(ierr);
471883469d8SStefano Zampini       local_size += subset_size;
472883469d8SStefano Zampini     }
473883469d8SStefano Zampini     ierr = PetscFree(dummy_idx);CHKERRQ(ierr);
4745db18549SStefano Zampini   }
4755db18549SStefano Zampini   ierr = MatAssemblyBegin(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
4765db18549SStefano Zampini   ierr = MatAssemblyEnd(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
4775db18549SStefano Zampini 
4785db18549SStefano Zampini   /* Global matrix of all assembled Schur on subsets */
479883469d8SStefano 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);
4805db18549SStefano Zampini   ierr = ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)sub_schurs->l2gmap),1,local_size,all_local_idx_G,PETSC_COPY_VALUES,&l2gmap_subsets);CHKERRQ(ierr);
4815db18549SStefano Zampini   ierr = MatCreateIS(PetscObjectComm((PetscObject)sub_schurs->l2gmap),1,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size,l2gmap_subsets,&work_mat);CHKERRQ(ierr);
4825db18549SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&l2gmap_subsets);CHKERRQ(ierr);
4835db18549SStefano Zampini   ierr = MatISSetLocalMat(work_mat,sub_schurs->S_Ej_all);CHKERRQ(ierr);
4845db18549SStefano Zampini   ierr = MatISGetMPIXAIJ(work_mat,MAT_INITIAL_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
4855db18549SStefano Zampini 
4865db18549SStefano Zampini   /* Get local part of (\sum_j S_Ej) */
487883469d8SStefano Zampini   ierr = PetscFree(all_local_idx_N);CHKERRQ(ierr);
488883469d8SStefano Zampini   ierr = ISCreateGeneral(PetscObjectComm((PetscObject)sub_schurs->l2gmap),local_size,all_local_idx_G,PETSC_OWN_POINTER,&temp_is);CHKERRQ(ierr);
489883469d8SStefano Zampini   ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
4905db18549SStefano Zampini   ierr = MatDestroy(&global_schur_subsets);CHKERRQ(ierr);
4913202ece2SStefano Zampini 
4923202ece2SStefano Zampini   /* Computation of Stildas */
4933202ece2SStefano Zampini   if (S_all) {
4943202ece2SStefano Zampini     PetscInt    n_all;
4953202ece2SStefano Zampini     PetscScalar *vals;
4963202ece2SStefano Zampini 
4973202ece2SStefano Zampini     /* Work arrays */
4983202ece2SStefano Zampini     ierr = MatGetSize(S_all,&n_all,NULL);CHKERRQ(ierr);
4993202ece2SStefano Zampini     ierr = MatDenseGetArray(S_all,&vals);CHKERRQ(ierr);
5003202ece2SStefano Zampini     ierr = MatDenseRestoreArray(S_all,&vals);CHKERRQ(ierr);
5013202ece2SStefano Zampini   }
5023202ece2SStefano Zampini 
5033202ece2SStefano Zampini   ierr = MatDestroy(&S_all);CHKERRQ(ierr);
5043202ece2SStefano Zampini   ierr = MatDestroy(&work_mat);CHKERRQ(ierr);
5055db18549SStefano Zampini   ierr = ISDestroy(&temp_is);CHKERRQ(ierr);
506b1b3d7a2SStefano Zampini   PetscFunctionReturn(0);
507b1b3d7a2SStefano Zampini }
508b1b3d7a2SStefano Zampini 
509b1b3d7a2SStefano Zampini #undef __FUNCT__
510b1b3d7a2SStefano Zampini #define __FUNCT__ "PCBDDCSubSchursInit"
5115db18549SStefano Zampini PetscErrorCode PCBDDCSubSchursInit(PCBDDCSubSchurs sub_schurs, Mat A, Mat S, IS is_I, IS is_B, PCBDDCGraph graph, ISLocalToGlobalMapping BtoNmap, PetscInt seqthreshold)
512b1b3d7a2SStefano Zampini {
513b1b3d7a2SStefano Zampini   IS                  *faces,*edges,*all_cc;
514b1b3d7a2SStefano Zampini   PetscInt            *index_sequential,*index_parallel;
515b1b3d7a2SStefano Zampini   PetscInt            *auxlocal_sequential,*auxlocal_parallel;
516b1b3d7a2SStefano Zampini   PetscInt            *auxglobal_sequential,*auxglobal_parallel;
517b1b3d7a2SStefano Zampini   PetscInt            *auxmapping;//,*idxs;
518b1b3d7a2SStefano Zampini   PetscInt            i,max_subset_size;
519b1b3d7a2SStefano Zampini   PetscInt            n_sequential_problems,n_local_sequential_problems,n_parallel_problems,n_local_parallel_problems;
520b1b3d7a2SStefano Zampini   PetscInt            n_faces,n_edges,n_all_cc;
521b1b3d7a2SStefano Zampini   PetscBool           is_sorted;
522b1b3d7a2SStefano Zampini   PetscErrorCode  ierr;
523b1b3d7a2SStefano Zampini 
524b1b3d7a2SStefano Zampini   PetscFunctionBegin;
525b1b3d7a2SStefano Zampini   ierr = ISSorted(is_I,&is_sorted);CHKERRQ(ierr);
526b1b3d7a2SStefano Zampini   if (!is_sorted) {
527b1b3d7a2SStefano Zampini     SETERRQ(PetscObjectComm((PetscObject)is_I),PETSC_ERR_PLIB,"IS for I dofs should be shorted");
528b1b3d7a2SStefano Zampini   }
529b1b3d7a2SStefano Zampini   ierr = ISSorted(is_B,&is_sorted);CHKERRQ(ierr);
530b1b3d7a2SStefano Zampini   if (!is_sorted) {
531b1b3d7a2SStefano Zampini     SETERRQ(PetscObjectComm((PetscObject)is_B),PETSC_ERR_PLIB,"IS for B dofs should be shorted");
532b1b3d7a2SStefano Zampini   }
533b1b3d7a2SStefano Zampini 
534b1b3d7a2SStefano Zampini   /* reset any previous data */
535b1b3d7a2SStefano Zampini   ierr = PCBDDCSubSchursReset(sub_schurs);CHKERRQ(ierr);
536b1b3d7a2SStefano Zampini 
537b1b3d7a2SStefano Zampini   /* get index sets for faces and edges */
538b1b3d7a2SStefano Zampini   ierr = PCBDDCGraphGetCandidatesIS(graph,&n_faces,&faces,&n_edges,&edges,NULL);CHKERRQ(ierr);
539b1b3d7a2SStefano Zampini   n_all_cc = n_faces+n_edges;
540b1b3d7a2SStefano Zampini   ierr = PetscMalloc1(n_all_cc,&all_cc);CHKERRQ(ierr);
541b1b3d7a2SStefano Zampini   for (i=0;i<n_faces;i++) {
542b1b3d7a2SStefano Zampini     all_cc[i] = faces[i];
543b1b3d7a2SStefano Zampini   }
544b1b3d7a2SStefano Zampini   for (i=0;i<n_edges;i++) {
545b1b3d7a2SStefano Zampini     all_cc[n_faces+i] = edges[i];
546b1b3d7a2SStefano Zampini   }
547b1b3d7a2SStefano Zampini   ierr = PetscFree(faces);CHKERRQ(ierr);
548b1b3d7a2SStefano Zampini   ierr = PetscFree(edges);CHKERRQ(ierr);
549b1b3d7a2SStefano Zampini 
550b1b3d7a2SStefano Zampini   /* map interface's subsets */
551b1b3d7a2SStefano Zampini   max_subset_size = 0;
552b1b3d7a2SStefano Zampini   for (i=0;i<n_all_cc;i++) {
553b1b3d7a2SStefano Zampini     PetscInt subset_size;
554b1b3d7a2SStefano Zampini     ierr = ISGetLocalSize(all_cc[i],&subset_size);CHKERRQ(ierr);
555b1b3d7a2SStefano Zampini     max_subset_size = PetscMax(max_subset_size,subset_size);
556b1b3d7a2SStefano Zampini   }
557b1b3d7a2SStefano Zampini   ierr = PetscMalloc1(max_subset_size,&auxmapping);CHKERRQ(ierr);
558b1b3d7a2SStefano Zampini   ierr = PetscMalloc2(graph->ncc,&auxlocal_sequential,
559b1b3d7a2SStefano Zampini                       graph->ncc,&auxlocal_parallel);CHKERRQ(ierr);
560b1b3d7a2SStefano Zampini   ierr = PetscMalloc2(graph->ncc,&index_sequential,
561b1b3d7a2SStefano Zampini                       graph->ncc,&index_parallel);CHKERRQ(ierr);
562b1b3d7a2SStefano Zampini 
563883469d8SStefano Zampini   /* if threshold is negative uses all sequential problems (possibly using MUMPS) */
564883469d8SStefano Zampini   sub_schurs->use_mumps = PETSC_FALSE;
565883469d8SStefano Zampini   if (seqthreshold < 0) {
566883469d8SStefano Zampini     seqthreshold = max_subset_size;
567883469d8SStefano Zampini #if defined(PETSC_HAVE_MUMPS)
568883469d8SStefano Zampini     sub_schurs->use_mumps = !!A;
569883469d8SStefano Zampini #endif
570883469d8SStefano Zampini   }
571883469d8SStefano Zampini 
572b1b3d7a2SStefano Zampini 
573b1b3d7a2SStefano Zampini   /* determine which problem has to be solved in parallel or sequentially */
574b1b3d7a2SStefano Zampini   n_local_sequential_problems = 0;
575b1b3d7a2SStefano Zampini   n_local_parallel_problems = 0;
576b1b3d7a2SStefano Zampini   for (i=0;i<n_all_cc;i++) {
577b1b3d7a2SStefano Zampini     PetscInt       subset_size,j,min_loc = 0;
578b1b3d7a2SStefano Zampini     const PetscInt *idxs;
579b1b3d7a2SStefano Zampini 
580b1b3d7a2SStefano Zampini     ierr = ISGetLocalSize(all_cc[i],&subset_size);CHKERRQ(ierr);
581b1b3d7a2SStefano Zampini     ierr = ISGetIndices(all_cc[i],&idxs);CHKERRQ(ierr);
582b1b3d7a2SStefano Zampini     ierr = ISLocalToGlobalMappingApply(graph->l2gmap,subset_size,idxs,auxmapping);CHKERRQ(ierr);
583b1b3d7a2SStefano Zampini     for (j=1;j<subset_size;j++) {
584b1b3d7a2SStefano Zampini       if (auxmapping[j]<auxmapping[min_loc]) {
585b1b3d7a2SStefano Zampini         min_loc = j;
586b1b3d7a2SStefano Zampini       }
587b1b3d7a2SStefano Zampini     }
588b1b3d7a2SStefano Zampini     if (subset_size > seqthreshold) {
589b1b3d7a2SStefano Zampini       index_parallel[n_local_parallel_problems] = i;
590b1b3d7a2SStefano Zampini       auxlocal_parallel[n_local_parallel_problems] = idxs[min_loc];
591b1b3d7a2SStefano Zampini       n_local_parallel_problems++;
592b1b3d7a2SStefano Zampini     } else {
593b1b3d7a2SStefano Zampini       index_sequential[n_local_sequential_problems] = i;
594b1b3d7a2SStefano Zampini       auxlocal_sequential[n_local_sequential_problems] = idxs[min_loc];
595b1b3d7a2SStefano Zampini       n_local_sequential_problems++;
596b1b3d7a2SStefano Zampini     }
597b1b3d7a2SStefano Zampini     ierr = ISRestoreIndices(all_cc[i],&idxs);CHKERRQ(ierr);
598b1b3d7a2SStefano Zampini   }
599b1b3d7a2SStefano Zampini 
600b1b3d7a2SStefano Zampini   /* Number parallel problems */
601b1b3d7a2SStefano Zampini   auxglobal_parallel = 0;
602b1b3d7a2SStefano Zampini   ierr = PCBDDCSubsetNumbering(PetscObjectComm((PetscObject)graph->l2gmap),graph->l2gmap,n_local_parallel_problems,auxlocal_parallel,PETSC_NULL,&n_parallel_problems,&auxglobal_parallel);CHKERRQ(ierr);
603b1b3d7a2SStefano Zampini 
604b1b3d7a2SStefano Zampini   /* Number sequential problems */
605b1b3d7a2SStefano Zampini   auxglobal_sequential = 0;
606b1b3d7a2SStefano Zampini   ierr = PCBDDCSubsetNumbering(PetscObjectComm((PetscObject)graph->l2gmap),graph->l2gmap,n_local_sequential_problems,auxlocal_sequential,PETSC_NULL,&n_sequential_problems,&auxglobal_sequential);CHKERRQ(ierr);
607b1b3d7a2SStefano Zampini 
608b1b3d7a2SStefano Zampini   /* update info in sub_schurs */
609883469d8SStefano Zampini   if (sub_schurs->use_mumps && A) {
6101e9c79c2SStefano Zampini     ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
6111e9c79c2SStefano Zampini     sub_schurs->A = A;
6121e9c79c2SStefano Zampini   }
613b1b3d7a2SStefano Zampini   ierr = PetscObjectReference((PetscObject)S);CHKERRQ(ierr);
614b1b3d7a2SStefano Zampini   sub_schurs->S = S;
615b1b3d7a2SStefano Zampini   ierr = PetscObjectReference((PetscObject)is_I);CHKERRQ(ierr);
616b1b3d7a2SStefano Zampini   sub_schurs->is_I = is_I;
617b1b3d7a2SStefano Zampini   ierr = PetscObjectReference((PetscObject)is_B);CHKERRQ(ierr);
618b1b3d7a2SStefano Zampini   sub_schurs->is_B = is_B;
6195db18549SStefano Zampini   ierr = PetscObjectReference((PetscObject)graph->l2gmap);CHKERRQ(ierr);
6205db18549SStefano Zampini   sub_schurs->l2gmap = graph->l2gmap;
6215db18549SStefano Zampini   ierr = PetscObjectReference((PetscObject)BtoNmap);CHKERRQ(ierr);
6225db18549SStefano Zampini   sub_schurs->BtoNmap = BtoNmap;
623b1b3d7a2SStefano Zampini   sub_schurs->n_subs_seq = n_local_sequential_problems;
624b1b3d7a2SStefano Zampini   sub_schurs->n_subs_par = n_local_parallel_problems;
625b1b3d7a2SStefano Zampini   sub_schurs->n_subs_seq_g = n_sequential_problems;
626b1b3d7a2SStefano Zampini   sub_schurs->n_subs_par_g = n_parallel_problems;
627b1b3d7a2SStefano Zampini   sub_schurs->n_subs = sub_schurs->n_subs_seq + sub_schurs->n_subs_par;
628b1b3d7a2SStefano Zampini   sub_schurs->is_subs = all_cc;
629b1b3d7a2SStefano Zampini   sub_schurs->index_sequential = index_sequential;
630b1b3d7a2SStefano Zampini   sub_schurs->index_parallel = index_parallel;
631b1b3d7a2SStefano Zampini   sub_schurs->auxglobal_sequential = auxglobal_sequential;
632b1b3d7a2SStefano Zampini   sub_schurs->auxglobal_parallel = auxglobal_parallel;
633b1b3d7a2SStefano Zampini 
634b1b3d7a2SStefano Zampini   /* free workspace */
635b1b3d7a2SStefano Zampini   ierr = PetscFree(auxmapping);CHKERRQ(ierr);
636b1b3d7a2SStefano Zampini   ierr = PetscFree2(auxlocal_sequential,auxlocal_parallel);CHKERRQ(ierr);
637b1b3d7a2SStefano Zampini 
638b1b3d7a2SStefano Zampini   PetscFunctionReturn(0);
639b1b3d7a2SStefano Zampini }
640b1b3d7a2SStefano Zampini 
641b1b3d7a2SStefano Zampini #undef __FUNCT__
64234a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCSubSchursCreate"
64334a97f8cSStefano Zampini PetscErrorCode PCBDDCSubSchursCreate(PCBDDCSubSchurs *sub_schurs)
64434a97f8cSStefano Zampini {
64534a97f8cSStefano Zampini   PCBDDCSubSchurs schurs_ctx;
64634a97f8cSStefano Zampini   PetscErrorCode  ierr;
64734a97f8cSStefano Zampini 
64834a97f8cSStefano Zampini   PetscFunctionBegin;
64934a97f8cSStefano Zampini   ierr = PetscNew(&schurs_ctx);CHKERRQ(ierr);
6505ff63025SStefano Zampini   schurs_ctx->n_subs = 0;
65134a97f8cSStefano Zampini   *sub_schurs = schurs_ctx;
65234a97f8cSStefano Zampini   PetscFunctionReturn(0);
65334a97f8cSStefano Zampini }
65434a97f8cSStefano Zampini 
65534a97f8cSStefano Zampini #undef __FUNCT__
65634a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCSubSchursDestroy"
65734a97f8cSStefano Zampini PetscErrorCode PCBDDCSubSchursDestroy(PCBDDCSubSchurs *sub_schurs)
65834a97f8cSStefano Zampini {
65934a97f8cSStefano Zampini   PetscErrorCode ierr;
66034a97f8cSStefano Zampini 
66134a97f8cSStefano Zampini   PetscFunctionBegin;
66234a97f8cSStefano Zampini   ierr = PCBDDCSubSchursReset(*sub_schurs);CHKERRQ(ierr);
66334a97f8cSStefano Zampini   ierr = PetscFree(*sub_schurs);CHKERRQ(ierr);
66434a97f8cSStefano Zampini   PetscFunctionReturn(0);
66534a97f8cSStefano Zampini }
66634a97f8cSStefano Zampini 
66734a97f8cSStefano Zampini #undef __FUNCT__
66834a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCSubSchursReset"
66934a97f8cSStefano Zampini PetscErrorCode PCBDDCSubSchursReset(PCBDDCSubSchurs sub_schurs)
67034a97f8cSStefano Zampini {
67134a97f8cSStefano Zampini   PetscInt       i;
67234a97f8cSStefano Zampini   PetscErrorCode ierr;
67334a97f8cSStefano Zampini 
67434a97f8cSStefano Zampini   PetscFunctionBegin;
6751e9c79c2SStefano Zampini   ierr = MatDestroy(&sub_schurs->A);CHKERRQ(ierr);
676b1b3d7a2SStefano Zampini   ierr = MatDestroy(&sub_schurs->S);CHKERRQ(ierr);
677b1b3d7a2SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_I);CHKERRQ(ierr);
678b1b3d7a2SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_B);CHKERRQ(ierr);
6795db18549SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->l2gmap);CHKERRQ(ierr);
6805db18549SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->BtoNmap);CHKERRQ(ierr);
68141c3ba1bSStefano Zampini   ierr = MatDestroy(&sub_schurs->S_Ej_all);CHKERRQ(ierr);
68241c3ba1bSStefano Zampini   ierr = MatDestroy(&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
6835db18549SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_Ej_all);CHKERRQ(ierr);
68434a97f8cSStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
685b1b3d7a2SStefano Zampini     ierr = ISDestroy(&sub_schurs->is_subs[i]);CHKERRQ(ierr);
68634a97f8cSStefano Zampini     ierr = ISDestroy(&sub_schurs->is_AEj_B[i]);CHKERRQ(ierr);
68734a97f8cSStefano Zampini     ierr = MatDestroy(&sub_schurs->S_Ej[i]);CHKERRQ(ierr);
6889221af80SStefano Zampini     ierr = VecDestroy(&sub_schurs->work1[i]);CHKERRQ(ierr);
6899221af80SStefano Zampini     ierr = VecDestroy(&sub_schurs->work2[i]);CHKERRQ(ierr);
69034a97f8cSStefano Zampini   }
69168270318SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_I_layer);CHKERRQ(ierr);
6925ff63025SStefano Zampini   if (sub_schurs->n_subs) {
693b1b3d7a2SStefano Zampini     ierr = PetscFree(sub_schurs->is_subs);CHKERRQ(ierr);
69468270318SStefano Zampini     ierr = PetscFree4(sub_schurs->is_AEj_B,sub_schurs->S_Ej,sub_schurs->work1,sub_schurs->work2);CHKERRQ(ierr);
695b1b3d7a2SStefano Zampini     ierr = PetscFree2(sub_schurs->index_sequential,sub_schurs->index_parallel);CHKERRQ(ierr);
696b1b3d7a2SStefano Zampini     ierr = PetscFree(sub_schurs->auxglobal_sequential);CHKERRQ(ierr);
697b1b3d7a2SStefano Zampini     ierr = PetscFree(sub_schurs->auxglobal_parallel);CHKERRQ(ierr);
6985ff63025SStefano Zampini   }
69934a97f8cSStefano Zampini   sub_schurs->n_subs = 0;
70034a97f8cSStefano Zampini   PetscFunctionReturn(0);
70134a97f8cSStefano Zampini }
70234a97f8cSStefano Zampini 
70334a97f8cSStefano Zampini #undef __FUNCT__
70434a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCAdjGetNextLayer_Private"
7052a155e38SStefano Zampini PETSC_STATIC_INLINE PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt* queue_tip,PetscInt n_prev,PetscBT touched,PetscInt* xadj,PetscInt* adjncy,PetscInt* n_added)
70634a97f8cSStefano Zampini {
70734a97f8cSStefano Zampini   PetscInt       i,j,n;
70834a97f8cSStefano Zampini   PetscErrorCode ierr;
70934a97f8cSStefano Zampini 
71034a97f8cSStefano Zampini   PetscFunctionBegin;
71134a97f8cSStefano Zampini   n = 0;
71234a97f8cSStefano Zampini   for (i=-n_prev;i<0;i++) {
71334a97f8cSStefano Zampini     PetscInt start_dof = queue_tip[i];
71434a97f8cSStefano Zampini     for (j=xadj[start_dof];j<xadj[start_dof+1];j++) {
71534a97f8cSStefano Zampini       PetscInt dof = adjncy[j];
71634a97f8cSStefano Zampini       if (!PetscBTLookup(touched,dof)) {
71734a97f8cSStefano Zampini         ierr = PetscBTSet(touched,dof);CHKERRQ(ierr);
71834a97f8cSStefano Zampini         queue_tip[n] = dof;
71934a97f8cSStefano Zampini         n++;
72034a97f8cSStefano Zampini       }
72134a97f8cSStefano Zampini     }
72234a97f8cSStefano Zampini   }
72334a97f8cSStefano Zampini   *n_added = n;
72434a97f8cSStefano Zampini   PetscFunctionReturn(0);
72534a97f8cSStefano Zampini }
726