xref: /petsc/src/ksp/pc/impls/bddc/bddcschurs.c (revision 3927de2ed0f0aec79c45dc8add25bfe718efd5fa)
134a97f8cSStefano Zampini #include <../src/ksp/pc/impls/bddc/bddc.h>
234a97f8cSStefano Zampini #include <../src/ksp/pc/impls/bddc/bddcprivate.h>
308122e43SStefano Zampini #include <petscblaslapack.h>
434a97f8cSStefano Zampini 
53202ece2SStefano Zampini PETSC_STATIC_INLINE PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt*,PetscInt,PetscBT,PetscInt*,PetscInt*,PetscInt*);
63202ece2SStefano Zampini static PetscErrorCode PCBDDCComputeExplicitSchur(Mat M, Mat *S);
73202ece2SStefano Zampini 
83202ece2SStefano Zampini #undef __FUNCT__
93202ece2SStefano Zampini #define __FUNCT__ "PCBDDCComputeExplicitSchur"
103202ece2SStefano Zampini static PetscErrorCode PCBDDCComputeExplicitSchur(Mat M, Mat *S)
113202ece2SStefano Zampini {
123202ece2SStefano Zampini   Mat            B, C, D, Bd, Cd, AinvBd;
133202ece2SStefano Zampini   KSP            ksp;
143202ece2SStefano Zampini   PC             pc;
153202ece2SStefano Zampini   PetscBool      isLU, isILU, isCHOL, Bdense, Cdense;
163202ece2SStefano Zampini   PetscReal      fill = 2.0;
173202ece2SStefano Zampini   PetscMPIInt    size;
183202ece2SStefano Zampini   PetscErrorCode ierr;
193202ece2SStefano Zampini 
203202ece2SStefano Zampini   PetscFunctionBegin;
213202ece2SStefano Zampini   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)M),&size);CHKERRQ(ierr);
223202ece2SStefano Zampini   if (size != 1) {
233202ece2SStefano Zampini     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not for parallel matrices");
243202ece2SStefano Zampini   }
253202ece2SStefano Zampini   ierr = MatSchurComplementGetSubMatrices(M, NULL, NULL, &B, &C, &D);CHKERRQ(ierr);
263202ece2SStefano Zampini   ierr = MatSchurComplementGetKSP(M, &ksp);CHKERRQ(ierr);
273202ece2SStefano Zampini   ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr);
283202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) pc, PCLU, &isLU);CHKERRQ(ierr);
293202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) pc, PCILU, &isILU);CHKERRQ(ierr);
303202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) pc, PCCHOLESKY, &isCHOL);CHKERRQ(ierr);
313202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) B, MATSEQDENSE, &Bdense);CHKERRQ(ierr);
323202ece2SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject) C, MATSEQDENSE, &Cdense);CHKERRQ(ierr);
333202ece2SStefano Zampini   if (!Bdense) {
343202ece2SStefano Zampini     ierr = MatConvert(B, MATSEQDENSE, MAT_INITIAL_MATRIX, &Bd);CHKERRQ(ierr);
353202ece2SStefano Zampini   } else {
363202ece2SStefano Zampini     Bd = B;
373202ece2SStefano Zampini   }
383202ece2SStefano Zampini 
393202ece2SStefano Zampini   if (isLU || isILU || isCHOL) {
403202ece2SStefano Zampini     Mat fact;
413202ece2SStefano Zampini 
423202ece2SStefano Zampini     ierr = KSPSetUp(ksp);CHKERRQ(ierr);
433202ece2SStefano Zampini     ierr = PCFactorGetMatrix(pc, &fact);CHKERRQ(ierr);
443202ece2SStefano Zampini     ierr = MatDuplicate(Bd, MAT_DO_NOT_COPY_VALUES, &AinvBd);CHKERRQ(ierr);
453202ece2SStefano Zampini     ierr = MatMatSolve(fact, Bd, AinvBd);CHKERRQ(ierr);
463202ece2SStefano Zampini   } else {
473202ece2SStefano Zampini     Mat Ainvd;
483202ece2SStefano Zampini 
493202ece2SStefano Zampini     ierr = PCComputeExplicitOperator(pc, &Ainvd);CHKERRQ(ierr);
503202ece2SStefano Zampini     ierr = MatMatMult(Ainvd, Bd, MAT_INITIAL_MATRIX, fill, &AinvBd);CHKERRQ(ierr);
513202ece2SStefano Zampini     ierr = MatDestroy(&Ainvd);CHKERRQ(ierr);
523202ece2SStefano Zampini   }
533202ece2SStefano Zampini   if (!Bdense) {
543202ece2SStefano Zampini     ierr = MatDestroy(&Bd);CHKERRQ(ierr);
553202ece2SStefano Zampini   }
563202ece2SStefano Zampini   if (!Cdense) {
573202ece2SStefano Zampini     ierr = MatConvert(C, MATSEQDENSE, MAT_INITIAL_MATRIX, &Cd);CHKERRQ(ierr);
583202ece2SStefano Zampini   } else {
593202ece2SStefano Zampini     Cd = C;
603202ece2SStefano Zampini   }
613202ece2SStefano Zampini 
623202ece2SStefano Zampini   ierr = MatMatMult(Cd, AinvBd, MAT_INITIAL_MATRIX, fill, S);CHKERRQ(ierr);
633202ece2SStefano Zampini   ierr = MatDestroy(&AinvBd);CHKERRQ(ierr);
643202ece2SStefano Zampini   if (!Cdense) {
653202ece2SStefano Zampini     ierr = MatDestroy(&Cd);CHKERRQ(ierr);
663202ece2SStefano Zampini   }
673202ece2SStefano Zampini 
683202ece2SStefano Zampini   if (D) {
693202ece2SStefano Zampini     Mat       Dd;
703202ece2SStefano Zampini     PetscBool Ddense;
713202ece2SStefano Zampini 
723202ece2SStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)D,MATSEQDENSE,&Ddense);CHKERRQ(ierr);
733202ece2SStefano Zampini     if (!Ddense) {
743202ece2SStefano Zampini       ierr = MatConvert(D, MATSEQDENSE, MAT_INITIAL_MATRIX, &Dd);CHKERRQ(ierr);
753202ece2SStefano Zampini     } else {
763202ece2SStefano Zampini       Dd = D;
773202ece2SStefano Zampini     }
783202ece2SStefano Zampini     ierr = MatAYPX(*S,-1.0,Dd,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
793202ece2SStefano Zampini     if (!Ddense) {
803202ece2SStefano Zampini       ierr = MatDestroy(&Dd);CHKERRQ(ierr);
813202ece2SStefano Zampini     }
823202ece2SStefano Zampini   } else {
833202ece2SStefano Zampini     ierr = MatScale(*S,-1.0);CHKERRQ(ierr);
843202ece2SStefano Zampini   }
853202ece2SStefano Zampini   PetscFunctionReturn(0);
863202ece2SStefano Zampini }
8734a97f8cSStefano Zampini 
8834a97f8cSStefano Zampini #undef __FUNCT__
891580ed26SStefano Zampini #define __FUNCT__ "PCBDDCSubSchursSetUp"
9008122e43SStefano Zampini PetscErrorCode PCBDDCSubSchursSetUp(PCBDDCSubSchurs sub_schurs, PetscInt xadj[], PetscInt adjncy[], PetscInt nlayers, PetscBool compute_Stilda, PetscBool deluxe, PetscBool use_edges, PetscBool use_faces)
91b1b3d7a2SStefano Zampini {
92b1b3d7a2SStefano Zampini   Mat                    A_II,A_IB,A_BI,A_BB;
93b1b3d7a2SStefano Zampini   Mat                    AE_II,*AE_IE,*AE_EI,*AE_EE;
94883469d8SStefano Zampini   Mat                    S_all,global_schur_subsets,work_mat;
9508122e43SStefano Zampini   Mat                    S_Ej_tilda_all,S_Ej_inv_all;
965db18549SStefano Zampini   ISLocalToGlobalMapping l2gmap_subsets;
975db18549SStefano Zampini   IS                     is_I,*is_subset_B,temp_is;
98883469d8SStefano Zampini   PetscInt               *nnz,*all_local_idx_G,*all_local_idx_B,*all_local_idx_N;
995db18549SStefano Zampini   PetscInt               i,subset_size,max_subset_size;
100883469d8SStefano Zampini   PetscInt               extra,local_size,global_size;
10108122e43SStefano Zampini   PetscBool              is_symmetric,Stilda_computed;
10208122e43SStefano Zampini   PetscBLASInt           B_N,B_ierr,B_lwork,*pivots;
10308122e43SStefano Zampini   PetscScalar            *work;
104b1b3d7a2SStefano Zampini   PetscErrorCode         ierr;
105b1b3d7a2SStefano Zampini 
106b1b3d7a2SStefano Zampini   PetscFunctionBegin;
107b1b3d7a2SStefano Zampini   /* get Schur complement matrices */
108883469d8SStefano Zampini   if (!sub_schurs->use_mumps) {
109b7eb3628SStefano Zampini     if (compute_Stilda) {
1109bb4a8caSStefano Zampini       SETERRQ(PetscObjectComm((PetscObject)sub_schurs->l2gmap),PETSC_ERR_SUP,"Computation of Stildas requires MUMPS");
111b7eb3628SStefano Zampini     }
112b1b3d7a2SStefano Zampini     ierr = MatSchurComplementGetSubMatrices(sub_schurs->S,&A_II,NULL,&A_IB,&A_BI,&A_BB);CHKERRQ(ierr);
113b1b3d7a2SStefano Zampini     ierr = PetscMalloc4(sub_schurs->n_subs,&is_subset_B,
114b1b3d7a2SStefano Zampini                         sub_schurs->n_subs,&AE_IE,
115b1b3d7a2SStefano Zampini                         sub_schurs->n_subs,&AE_EI,
116b1b3d7a2SStefano Zampini                         sub_schurs->n_subs,&AE_EE);CHKERRQ(ierr);
117b1b3d7a2SStefano Zampini   }
118b1b3d7a2SStefano Zampini 
119b1b3d7a2SStefano Zampini   /* determine interior problems */
120b96c3477SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_I_layer);CHKERRQ(ierr);
121b1b3d7a2SStefano Zampini   if (nlayers >= 0 && xadj != NULL && adjncy != NULL) { /* Interior problems can be different from the original one */
122b1b3d7a2SStefano Zampini     PetscBT                touched;
123b1b3d7a2SStefano Zampini     const PetscInt*        idx_B;
124b1b3d7a2SStefano Zampini     PetscInt               n_I,n_B,n_local_dofs,n_prev_added,j,layer,*local_numbering;
125b1b3d7a2SStefano Zampini 
126b1b3d7a2SStefano Zampini     /* get sizes */
127b1b3d7a2SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr);
128b1b3d7a2SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_B,&n_B);CHKERRQ(ierr);
129b1b3d7a2SStefano Zampini 
130b1b3d7a2SStefano Zampini     ierr = PetscMalloc1(n_I+n_B,&local_numbering);CHKERRQ(ierr);
131b1b3d7a2SStefano Zampini     ierr = PetscBTCreate(n_I+n_B,&touched);CHKERRQ(ierr);
132b1b3d7a2SStefano Zampini     ierr = PetscBTMemzero(n_I+n_B,touched);CHKERRQ(ierr);
133b1b3d7a2SStefano Zampini 
134b1b3d7a2SStefano Zampini     /* all boundary dofs must be skipped when adding layers */
135b1b3d7a2SStefano Zampini     ierr = ISGetIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr);
136b1b3d7a2SStefano Zampini     for (j=0;j<n_B;j++) {
137b1b3d7a2SStefano Zampini       ierr = PetscBTSet(touched,idx_B[j]);CHKERRQ(ierr);
138b1b3d7a2SStefano Zampini     }
139b1b3d7a2SStefano Zampini     ierr = PetscMemcpy(local_numbering,idx_B,n_B*sizeof(PetscInt));CHKERRQ(ierr);
140b1b3d7a2SStefano Zampini     ierr = ISRestoreIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr);
141b1b3d7a2SStefano Zampini 
142b1b3d7a2SStefano Zampini     /* add prescribed number of layers of dofs */
143b1b3d7a2SStefano Zampini     n_local_dofs = n_B;
144b1b3d7a2SStefano Zampini     n_prev_added = n_B;
145b1b3d7a2SStefano Zampini     for (layer=0;layer<nlayers;layer++) {
146b1b3d7a2SStefano Zampini       PetscInt n_added;
147b1b3d7a2SStefano Zampini       if (n_local_dofs == n_I+n_B) break;
148b1b3d7a2SStefano Zampini       if (n_local_dofs > n_I+n_B) {
149b1b3d7a2SStefano 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);
150b1b3d7a2SStefano Zampini       }
151b1b3d7a2SStefano Zampini       ierr = PCBDDCAdjGetNextLayer_Private(local_numbering+n_local_dofs,n_prev_added,touched,xadj,adjncy,&n_added);CHKERRQ(ierr);
152b1b3d7a2SStefano Zampini       n_prev_added = n_added;
153b1b3d7a2SStefano Zampini       n_local_dofs += n_added;
154b1b3d7a2SStefano Zampini       if (!n_added) break;
155b1b3d7a2SStefano Zampini     }
156b1b3d7a2SStefano Zampini     ierr = PetscBTDestroy(&touched);CHKERRQ(ierr);
157b1b3d7a2SStefano Zampini 
158883469d8SStefano Zampini     /* IS for I layer dofs in original numbering */
15968270318SStefano 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);
160b1b3d7a2SStefano Zampini     ierr = PetscFree(local_numbering);CHKERRQ(ierr);
16168270318SStefano Zampini     ierr = ISSort(sub_schurs->is_I_layer);CHKERRQ(ierr);
162883469d8SStefano Zampini     /* IS for I layer dofs in I numbering */
163883469d8SStefano Zampini     if (!sub_schurs->use_mumps) {
164b1b3d7a2SStefano Zampini       ISLocalToGlobalMapping ItoNmap;
165b1b3d7a2SStefano Zampini       ierr = ISLocalToGlobalMappingCreateIS(sub_schurs->is_I,&ItoNmap);CHKERRQ(ierr);
16668270318SStefano Zampini       ierr = ISGlobalToLocalMappingApplyIS(ItoNmap,IS_GTOLM_DROP,sub_schurs->is_I_layer,&is_I);CHKERRQ(ierr);
167b1b3d7a2SStefano Zampini       ierr = ISLocalToGlobalMappingDestroy(&ItoNmap);CHKERRQ(ierr);
168b1b3d7a2SStefano Zampini 
169b1b3d7a2SStefano Zampini       /* II block */
170b1b3d7a2SStefano Zampini       ierr = MatGetSubMatrix(A_II,is_I,is_I,MAT_INITIAL_MATRIX,&AE_II);CHKERRQ(ierr);
171b1b3d7a2SStefano Zampini     }
172b1b3d7a2SStefano Zampini   } else {
173b1b3d7a2SStefano Zampini     PetscInt n_I;
174b1b3d7a2SStefano Zampini 
175b1b3d7a2SStefano Zampini     /* IS for I dofs in original numbering */
176b1b3d7a2SStefano Zampini     ierr = PetscObjectReference((PetscObject)sub_schurs->is_I);CHKERRQ(ierr);
17768270318SStefano Zampini     sub_schurs->is_I_layer = sub_schurs->is_I;
178b1b3d7a2SStefano Zampini 
179b1b3d7a2SStefano Zampini     /* IS for I dofs in I numbering (strided 1) */
180883469d8SStefano Zampini     if (!sub_schurs->use_mumps) {
181b1b3d7a2SStefano Zampini       ierr = ISGetSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr);
182b1b3d7a2SStefano Zampini       ierr = ISCreateStride(PetscObjectComm((PetscObject)sub_schurs->is_I),n_I,0,1,&is_I);CHKERRQ(ierr);
183b1b3d7a2SStefano Zampini 
184b1b3d7a2SStefano Zampini       /* II block is the same */
185b1b3d7a2SStefano Zampini       ierr = PetscObjectReference((PetscObject)A_II);CHKERRQ(ierr);
186b1b3d7a2SStefano Zampini       AE_II = A_II;
187b1b3d7a2SStefano Zampini     }
188b1b3d7a2SStefano Zampini   }
189b1b3d7a2SStefano Zampini 
190883469d8SStefano Zampini   /* Get info on subset sizes and sum of all subsets sizes */
191883469d8SStefano Zampini   max_subset_size = 0;
192883469d8SStefano Zampini   local_size = 0;
193883469d8SStefano Zampini   for (i=0;i<sub_schurs->n_subs_seq;i++) {
194883469d8SStefano Zampini     PetscInt j = sub_schurs->index_sequential[i];
195b96c3477SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_subs[j],&subset_size);CHKERRQ(ierr);
196883469d8SStefano Zampini     max_subset_size = PetscMax(subset_size,max_subset_size);
197883469d8SStefano Zampini     local_size += subset_size;
198883469d8SStefano Zampini   }
199883469d8SStefano Zampini 
200883469d8SStefano Zampini   /* Work arrays for local indices */
201883469d8SStefano Zampini   ierr = PetscMalloc1(local_size,&all_local_idx_B);CHKERRQ(ierr);
202883469d8SStefano Zampini   extra = 0;
203883469d8SStefano Zampini   if (sub_schurs->use_mumps) {
204883469d8SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_I_layer,&extra);CHKERRQ(ierr);
205883469d8SStefano Zampini   }
206883469d8SStefano Zampini   ierr = PetscMalloc1(local_size+extra,&all_local_idx_N);CHKERRQ(ierr);
207883469d8SStefano Zampini   if (extra) {
208883469d8SStefano Zampini     const PetscInt *idxs;
209883469d8SStefano Zampini     ierr = ISGetIndices(sub_schurs->is_I_layer,&idxs);CHKERRQ(ierr);
210883469d8SStefano Zampini     ierr = PetscMemcpy(all_local_idx_N,idxs,extra*sizeof(PetscInt));CHKERRQ(ierr);
211883469d8SStefano Zampini     ierr = ISRestoreIndices(sub_schurs->is_I_layer,&idxs);CHKERRQ(ierr);
212883469d8SStefano Zampini   }
213883469d8SStefano Zampini   ierr = PetscMalloc1(local_size,&nnz);CHKERRQ(ierr);
214883469d8SStefano Zampini 
215883469d8SStefano Zampini   /* Get local indices in local numbering */
216883469d8SStefano Zampini   local_size = 0;
217883469d8SStefano Zampini   for (i=0;i<sub_schurs->n_subs_seq;i++) {
218883469d8SStefano Zampini     PetscInt j;
219883469d8SStefano Zampini     const    PetscInt *idxs;
220883469d8SStefano Zampini 
221883469d8SStefano Zampini     PetscInt local_problem_index = sub_schurs->index_sequential[i];
222b96c3477SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_subs[local_problem_index],&subset_size);CHKERRQ(ierr);
223b96c3477SStefano Zampini     ierr = ISGetIndices(sub_schurs->is_subs[local_problem_index],&idxs);CHKERRQ(ierr);
224883469d8SStefano Zampini     /* subset indices in local numbering */
225883469d8SStefano Zampini     ierr = PetscMemcpy(all_local_idx_N+local_size+extra,idxs,subset_size*sizeof(PetscInt));CHKERRQ(ierr);
226b96c3477SStefano Zampini     ierr = ISRestoreIndices(sub_schurs->is_subs[local_problem_index],&idxs);CHKERRQ(ierr);
227883469d8SStefano Zampini     for (j=0;j<subset_size;j++) nnz[local_size+j] = subset_size;
228883469d8SStefano Zampini     local_size += subset_size;
229883469d8SStefano Zampini   }
230883469d8SStefano Zampini 
231883469d8SStefano Zampini   S_all = NULL;
23208122e43SStefano Zampini   ierr = MatIsSymmetric(sub_schurs->A,0.0,&is_symmetric);CHKERRQ(ierr);
233883469d8SStefano Zampini   if (!sub_schurs->use_mumps) {
234883469d8SStefano Zampini     /* subsets in original and boundary numbering */
235883469d8SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
236b96c3477SStefano Zampini       ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,sub_schurs->is_subs[i],&is_subset_B[i]);CHKERRQ(ierr);
237b1b3d7a2SStefano Zampini     }
238b1b3d7a2SStefano Zampini 
239b1b3d7a2SStefano Zampini     /* EE block */
240b1b3d7a2SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
241b1b3d7a2SStefano Zampini       ierr = MatGetSubMatrix(A_BB,is_subset_B[i],is_subset_B[i],MAT_INITIAL_MATRIX,&AE_EE[i]);CHKERRQ(ierr);
242b1b3d7a2SStefano Zampini     }
243b1b3d7a2SStefano Zampini     /* IE block */
244b1b3d7a2SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
245b1b3d7a2SStefano Zampini       ierr = MatGetSubMatrix(A_IB,is_I,is_subset_B[i],MAT_INITIAL_MATRIX,&AE_IE[i]);CHKERRQ(ierr);
246b1b3d7a2SStefano Zampini     }
247b1b3d7a2SStefano Zampini     /* EI block */
248b1b3d7a2SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
249b1b3d7a2SStefano Zampini       ierr = MatGetSubMatrix(A_BI,is_subset_B[i],is_I,MAT_INITIAL_MATRIX,&AE_EI[i]);CHKERRQ(ierr);
250b1b3d7a2SStefano Zampini     }
251b1b3d7a2SStefano Zampini 
252b1b3d7a2SStefano Zampini     /* setup Schur complements on subset */
253b1b3d7a2SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
254b96c3477SStefano Zampini       ierr = MatDestroy(&sub_schurs->S_Ej[i]);CHKERRQ(ierr);
255b1b3d7a2SStefano Zampini       ierr = MatCreateSchurComplement(AE_II,AE_II,AE_IE[i],AE_EI[i],AE_EE[i],&sub_schurs->S_Ej[i]);CHKERRQ(ierr);
256b1b3d7a2SStefano Zampini       if (AE_II == A_II) { /* we can reuse the same ksp */
257b1b3d7a2SStefano Zampini         KSP ksp;
258b1b3d7a2SStefano Zampini         ierr = MatSchurComplementGetKSP(sub_schurs->S,&ksp);CHKERRQ(ierr);
259b1b3d7a2SStefano Zampini         ierr = MatSchurComplementSetKSP(sub_schurs->S_Ej[i],ksp);CHKERRQ(ierr);
260b1b3d7a2SStefano Zampini       } else { /* build new ksp object which inherits ksp and pc types from the original one */
261b1b3d7a2SStefano Zampini         KSP      origksp,schurksp;
262b1b3d7a2SStefano Zampini         PC       origpc,schurpc;
263b1b3d7a2SStefano Zampini         KSPType  ksp_type;
264b1b3d7a2SStefano Zampini         PCType   pc_type;
265b1b3d7a2SStefano Zampini         PetscInt n_internal;
266b1b3d7a2SStefano Zampini 
267b1b3d7a2SStefano Zampini         ierr = MatSchurComplementGetKSP(sub_schurs->S,&origksp);CHKERRQ(ierr);
268b1b3d7a2SStefano Zampini         ierr = MatSchurComplementGetKSP(sub_schurs->S_Ej[i],&schurksp);CHKERRQ(ierr);
269b1b3d7a2SStefano Zampini         ierr = KSPGetType(origksp,&ksp_type);CHKERRQ(ierr);
270b1b3d7a2SStefano Zampini         ierr = KSPSetType(schurksp,ksp_type);CHKERRQ(ierr);
271b1b3d7a2SStefano Zampini         ierr = KSPGetPC(schurksp,&schurpc);CHKERRQ(ierr);
272b1b3d7a2SStefano Zampini         ierr = KSPGetPC(origksp,&origpc);CHKERRQ(ierr);
273b1b3d7a2SStefano Zampini         ierr = PCGetType(origpc,&pc_type);CHKERRQ(ierr);
274b1b3d7a2SStefano Zampini         ierr = PCSetType(schurpc,pc_type);CHKERRQ(ierr);
275b1b3d7a2SStefano Zampini         ierr = ISGetSize(is_I,&n_internal);CHKERRQ(ierr);
276b1b3d7a2SStefano Zampini         if (n_internal) { /* UMFPACK gives error with 0 sized problems */
277b1b3d7a2SStefano Zampini           MatSolverPackage solver=NULL;
278b1b3d7a2SStefano Zampini           ierr = PCFactorGetMatSolverPackage(origpc,(const MatSolverPackage*)&solver);CHKERRQ(ierr);
279b1b3d7a2SStefano Zampini           if (solver) {
280b1b3d7a2SStefano Zampini             ierr = PCFactorSetMatSolverPackage(schurpc,solver);CHKERRQ(ierr);
281b1b3d7a2SStefano Zampini           }
282b1b3d7a2SStefano Zampini         }
283b1b3d7a2SStefano Zampini         ierr = KSPSetUp(schurksp);CHKERRQ(ierr);
284b1b3d7a2SStefano Zampini       }
285b1b3d7a2SStefano Zampini     }
286b1b3d7a2SStefano Zampini     /* free */
287b1b3d7a2SStefano Zampini     ierr = ISDestroy(&is_I);CHKERRQ(ierr);
288b1b3d7a2SStefano Zampini     ierr = MatDestroy(&AE_II);CHKERRQ(ierr);
289b1b3d7a2SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
290b1b3d7a2SStefano Zampini       ierr = MatDestroy(&AE_EE[i]);CHKERRQ(ierr);
291b1b3d7a2SStefano Zampini       ierr = MatDestroy(&AE_IE[i]);CHKERRQ(ierr);
292b1b3d7a2SStefano Zampini       ierr = MatDestroy(&AE_EI[i]);CHKERRQ(ierr);
293b1b3d7a2SStefano Zampini       ierr = ISDestroy(&is_subset_B[i]);CHKERRQ(ierr);
294b1b3d7a2SStefano Zampini     }
295b1b3d7a2SStefano Zampini     ierr = PetscFree4(is_subset_B,AE_IE,AE_EI,AE_EE);CHKERRQ(ierr);
296883469d8SStefano Zampini #if defined(PETSC_HAVE_MUMPS)
297883469d8SStefano Zampini   } else {
298883469d8SStefano Zampini     Mat           A,F;
299883469d8SStefano Zampini     IS            is_A_all;
300883469d8SStefano Zampini     PetscInt      *idxs_schur,n_I;
301883469d8SStefano Zampini 
302883469d8SStefano Zampini     /* get working mat */
303883469d8SStefano Zampini     ierr = ISGetLocalSize(sub_schurs->is_I_layer,&n_I);CHKERRQ(ierr);
304883469d8SStefano Zampini     ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size+n_I,all_local_idx_N,PETSC_COPY_VALUES,&is_A_all);CHKERRQ(ierr);
305b96c3477SStefano Zampini     ierr = MatGetSubMatrixUnsorted(sub_schurs->A,is_A_all,is_A_all,MAT_INITIAL_MATRIX,&A);CHKERRQ(ierr);
306883469d8SStefano Zampini     ierr = ISDestroy(&is_A_all);CHKERRQ(ierr);
307883469d8SStefano Zampini 
30808122e43SStefano Zampini     if (n_I) {
309883469d8SStefano Zampini       if (is_symmetric) {
310883469d8SStefano Zampini         ierr = MatGetFactor(A,MATSOLVERMUMPS,MAT_FACTOR_CHOLESKY,&F);CHKERRQ(ierr);
311883469d8SStefano Zampini       } else {
312883469d8SStefano Zampini         ierr = MatGetFactor(A,MATSOLVERMUMPS,MAT_FACTOR_LU,&F);CHKERRQ(ierr);
313883469d8SStefano Zampini       }
314883469d8SStefano Zampini 
315883469d8SStefano Zampini       /* subsets ordered last */
316883469d8SStefano Zampini       ierr = PetscMalloc1(local_size,&idxs_schur);CHKERRQ(ierr);
317883469d8SStefano Zampini       for (i=0;i<local_size;i++) {
318883469d8SStefano Zampini         idxs_schur[i] = n_I+i+1;
319883469d8SStefano Zampini       }
320883469d8SStefano Zampini       ierr = MatMumpsSetSchurIndices(F,local_size,idxs_schur);CHKERRQ(ierr);
321883469d8SStefano Zampini       ierr = PetscFree(idxs_schur);CHKERRQ(ierr);
322883469d8SStefano Zampini 
323883469d8SStefano Zampini       /* factorization step */
324883469d8SStefano Zampini       if (is_symmetric) {
325883469d8SStefano Zampini         ierr = MatCholeskyFactorSymbolic(F,A,NULL,NULL);CHKERRQ(ierr);
326883469d8SStefano Zampini         ierr = MatCholeskyFactorNumeric(F,A,NULL);CHKERRQ(ierr);
327883469d8SStefano Zampini       } else {
328883469d8SStefano Zampini         ierr = MatLUFactorSymbolic(F,A,NULL,NULL,NULL);CHKERRQ(ierr);
329883469d8SStefano Zampini         ierr = MatLUFactorNumeric(F,A,NULL);CHKERRQ(ierr);
330883469d8SStefano Zampini       }
331883469d8SStefano Zampini 
332883469d8SStefano Zampini       /* get explicit Schur Complement computed during numeric factorization */
333883469d8SStefano Zampini       ierr = MatMumpsGetSchurComplement(F,&S_all);CHKERRQ(ierr);
334883469d8SStefano Zampini       ierr = MatDestroy(&F);CHKERRQ(ierr);
33508122e43SStefano Zampini     } else {
33608122e43SStefano Zampini       ierr = MatConvert(A,MATSEQDENSE,MAT_INITIAL_MATRIX,&S_all);CHKERRQ(ierr);
33708122e43SStefano Zampini     }
33808122e43SStefano Zampini     ierr = MatDestroy(&A);CHKERRQ(ierr);
339883469d8SStefano Zampini #endif
340b1b3d7a2SStefano Zampini   }
3415db18549SStefano Zampini 
3425db18549SStefano Zampini   if (!sub_schurs->n_subs_seq_g) {
3435db18549SStefano Zampini     PetscFunctionReturn(0);
3445db18549SStefano Zampini   }
3455db18549SStefano Zampini 
3465db18549SStefano Zampini   /* subset indices in local boundary numbering */
347b96c3477SStefano Zampini   if (!sub_schurs->is_Ej_all) {
348883469d8SStefano Zampini     ierr = ISGlobalToLocalMappingApply(sub_schurs->BtoNmap,IS_GTOLM_DROP,local_size,all_local_idx_N+extra,&subset_size,all_local_idx_B);CHKERRQ(ierr);
3495db18549SStefano Zampini     if (subset_size != local_size) {
3505db18549SStefano Zampini       SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in sub_schurs serial (BtoNmap)! %d != %d\n",subset_size,local_size);
3515db18549SStefano Zampini     }
3525db18549SStefano Zampini     ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size,all_local_idx_B,PETSC_OWN_POINTER,&sub_schurs->is_Ej_all);CHKERRQ(ierr);
353b96c3477SStefano Zampini   }
3545db18549SStefano Zampini 
3553202ece2SStefano Zampini   /* Local matrix of all local Schur on subsets transposed */
356b96c3477SStefano Zampini   if (!sub_schurs->S_Ej_all) {
3575db18549SStefano Zampini     ierr = MatCreate(PETSC_COMM_SELF,&sub_schurs->S_Ej_all);CHKERRQ(ierr);
3585db18549SStefano Zampini     ierr = MatSetSizes(sub_schurs->S_Ej_all,PETSC_DECIDE,PETSC_DECIDE,local_size,local_size);CHKERRQ(ierr);
3595db18549SStefano Zampini     ierr = MatSetType(sub_schurs->S_Ej_all,MATAIJ);CHKERRQ(ierr);
3605db18549SStefano Zampini     ierr = MatSeqAIJSetPreallocation(sub_schurs->S_Ej_all,0,nnz);CHKERRQ(ierr);
361b96c3477SStefano Zampini   } else {
362b96c3477SStefano Zampini     ierr = MatZeroEntries(sub_schurs->S_Ej_all);CHKERRQ(ierr);
363b96c3477SStefano Zampini   }
364a1337663SStefano Zampini 
365a1337663SStefano Zampini   S_Ej_tilda_all = 0;
36608122e43SStefano Zampini   S_Ej_inv_all = 0;
36708122e43SStefano Zampini   if (sub_schurs->n_subs && deluxe) { /* workspace needed only for GETRI */
36808122e43SStefano Zampini     PetscScalar lwork;
36908122e43SStefano Zampini 
37008122e43SStefano Zampini     B_lwork = -1;
37108122e43SStefano Zampini     ierr = PetscBLASIntCast(max_subset_size,&B_N);CHKERRQ(ierr);
37208122e43SStefano Zampini     ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
37308122e43SStefano Zampini     PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,work,&B_N,pivots,&lwork,&B_lwork,&B_ierr));
37408122e43SStefano Zampini     ierr = PetscFPTrapPop();CHKERRQ(ierr);
37508122e43SStefano Zampini     if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to GETRI Lapack routine %d",(int)B_ierr);
37608122e43SStefano Zampini     ierr = PetscBLASIntCast((PetscInt)lwork,&B_lwork);CHKERRQ(ierr);
37708122e43SStefano Zampini     ierr = PetscMalloc2(B_lwork,&work,max_subset_size,&pivots);CHKERRQ(ierr);
37808122e43SStefano Zampini   } else {
37908122e43SStefano Zampini     work = NULL;
38008122e43SStefano Zampini     pivots = NULL;
38108122e43SStefano Zampini   }
38208122e43SStefano Zampini 
38308122e43SStefano Zampini   ierr = PetscBTMemzero(sub_schurs->n_subs,sub_schurs->computed_Stilda_subs);CHKERRQ(ierr);
384883469d8SStefano Zampini   if (!sub_schurs->use_mumps) {
3855db18549SStefano Zampini     PetscScalar *fill_vals;
3865db18549SStefano Zampini     PetscInt    *dummy_idx;
3875db18549SStefano Zampini 
3885db18549SStefano Zampini     /* Work arrays */
3893202ece2SStefano Zampini     ierr = PetscMalloc1(max_subset_size,&dummy_idx);CHKERRQ(ierr);
3903202ece2SStefano Zampini     /* Loop on local problems end compute Schur complements explicitly */
3915db18549SStefano Zampini     local_size = 0;
3925db18549SStefano Zampini     for (i=0;i<sub_schurs->n_subs_seq;i++) {
3933202ece2SStefano Zampini       Mat       S_Ej_expl;
3943202ece2SStefano Zampini       PetscInt  j,lpi = sub_schurs->index_sequential[i];
3953202ece2SStefano Zampini       PetscBool Sdense;
3965db18549SStefano Zampini 
3973202ece2SStefano Zampini       ierr = PCBDDCComputeExplicitSchur(sub_schurs->S_Ej[lpi],&S_Ej_expl);CHKERRQ(ierr);
398b96c3477SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_subs[lpi],&subset_size);CHKERRQ(ierr);
3993202ece2SStefano Zampini       ierr = PetscObjectTypeCompare((PetscObject)S_Ej_expl,MATSEQDENSE,&Sdense);CHKERRQ(ierr);
4003202ece2SStefano Zampini       if (Sdense) {
4015db18549SStefano Zampini         for (j=0;j<subset_size;j++) {
4025db18549SStefano Zampini           dummy_idx[j]=local_size+j;
4035db18549SStefano Zampini         }
4043202ece2SStefano Zampini         ierr = MatDenseGetArray(S_Ej_expl,&fill_vals);CHKERRQ(ierr);
4055db18549SStefano Zampini         ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,fill_vals,INSERT_VALUES);CHKERRQ(ierr);
4063202ece2SStefano Zampini         ierr = MatDenseRestoreArray(S_Ej_expl,&fill_vals);CHKERRQ(ierr);
4073202ece2SStefano Zampini       } else {
4083202ece2SStefano Zampini         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not yet implemented for sparse matrices");
4095db18549SStefano Zampini       }
4103202ece2SStefano Zampini       local_size += subset_size;
41165d8bf0aSStefano Zampini       ierr = MatDestroy(&sub_schurs->S_Ej[lpi]);CHKERRQ(ierr);
41265d8bf0aSStefano Zampini       sub_schurs->S_Ej[lpi] = S_Ej_expl;
4133202ece2SStefano Zampini     }
4142b973097SStefano Zampini     /* Stildas are not computed without mumps */
4152b973097SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
416b96c3477SStefano Zampini       ierr = MatDestroy(&sub_schurs->S_Ej_tilda[i]);CHKERRQ(ierr);
4172b973097SStefano Zampini     }
4183202ece2SStefano Zampini     ierr = PetscFree(dummy_idx);CHKERRQ(ierr);
419883469d8SStefano Zampini   } else {
420883469d8SStefano Zampini     PetscInt  *dummy_idx,n_all;
421883469d8SStefano Zampini 
422a1337663SStefano Zampini     if (compute_Stilda) {
423a1337663SStefano Zampini       ierr = MatCreate(PETSC_COMM_SELF,&S_Ej_tilda_all);CHKERRQ(ierr);
424a1337663SStefano Zampini       ierr = MatSetSizes(S_Ej_tilda_all,PETSC_DECIDE,PETSC_DECIDE,local_size,local_size);CHKERRQ(ierr);
425a1337663SStefano Zampini       ierr = MatSetType(S_Ej_tilda_all,MATAIJ);CHKERRQ(ierr);
426a1337663SStefano Zampini       ierr = MatSeqAIJSetPreallocation(S_Ej_tilda_all,0,nnz);CHKERRQ(ierr);
42708122e43SStefano Zampini       if (deluxe) {
42808122e43SStefano Zampini         ierr = MatCreate(PETSC_COMM_SELF,&S_Ej_inv_all);CHKERRQ(ierr);
42908122e43SStefano Zampini         ierr = MatSetSizes(S_Ej_inv_all,PETSC_DECIDE,PETSC_DECIDE,local_size,local_size);CHKERRQ(ierr);
43008122e43SStefano Zampini         ierr = MatSetType(S_Ej_inv_all,MATAIJ);CHKERRQ(ierr);
43108122e43SStefano Zampini         ierr = MatSeqAIJSetPreallocation(S_Ej_inv_all,0,nnz);CHKERRQ(ierr);
432a1337663SStefano Zampini       }
43308122e43SStefano Zampini     }
43408122e43SStefano Zampini 
43508122e43SStefano Zampini     ierr = MatGetSize(S_all,&n_all,NULL);CHKERRQ(ierr);
43608122e43SStefano Zampini     local_size = 0;
437883469d8SStefano Zampini     /* Work arrays */
438883469d8SStefano Zampini     ierr = PetscMalloc1(max_subset_size,&dummy_idx);CHKERRQ(ierr);
439883469d8SStefano Zampini 
44008122e43SStefano Zampini     Stilda_computed = PETSC_FALSE;
44165d8bf0aSStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
44265d8bf0aSStefano Zampini       IS  is_E;
44365d8bf0aSStefano Zampini       PetscScalar *vals;
44465d8bf0aSStefano Zampini       PetscInt j;
44565d8bf0aSStefano Zampini 
446b96c3477SStefano Zampini       ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
447a1337663SStefano Zampini       for (j=0;j<subset_size;j++) {
448a1337663SStefano Zampini         dummy_idx[j]=local_size+j;
449a1337663SStefano Zampini       }
45065d8bf0aSStefano Zampini       ierr = ISCreateStride(PetscObjectComm((PetscObject)sub_schurs->is_I),subset_size,local_size,1,&is_E);CHKERRQ(ierr);
451b96c3477SStefano Zampini       if (sub_schurs->S_Ej[i]) {
452b96c3477SStefano Zampini         ierr = MatGetSubMatrix(S_all,is_E,is_E,MAT_REUSE_MATRIX,&sub_schurs->S_Ej[i]);CHKERRQ(ierr);
453b96c3477SStefano Zampini       } else {
454a1337663SStefano Zampini         ierr = MatGetSubMatrix(S_all,is_E,is_E,MAT_INITIAL_MATRIX,&sub_schurs->S_Ej[i]);CHKERRQ(ierr);
455b96c3477SStefano Zampini       }
456a1337663SStefano Zampini       ierr = MatDenseGetArray(sub_schurs->S_Ej[i],&vals);CHKERRQ(ierr);
457a1337663SStefano Zampini       ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,vals,INSERT_VALUES);CHKERRQ(ierr);
458a1337663SStefano Zampini       ierr = MatDenseRestoreArray(sub_schurs->S_Ej[i],&vals);CHKERRQ(ierr);
459a1337663SStefano Zampini 
4606da91a30SStefano Zampini       if (compute_Stilda && ((PetscBTLookup(sub_schurs->is_edge,i) && use_edges) || use_faces)) {
46108122e43SStefano Zampini         ierr = PetscBTSet(sub_schurs->computed_Stilda_subs,i);CHKERRQ(ierr);
46208122e43SStefano Zampini         Stilda_computed = PETSC_TRUE;
46308122e43SStefano Zampini       }
46408122e43SStefano Zampini 
46508122e43SStefano Zampini       if (PetscBTLookup(sub_schurs->computed_Stilda_subs,i)) {
46665d8bf0aSStefano Zampini         Mat Stilda;
46765d8bf0aSStefano Zampini         if (sub_schurs->n_subs == 1) {
46865d8bf0aSStefano Zampini           ierr = PetscObjectReference((PetscObject)sub_schurs->S_Ej[i]);CHKERRQ(ierr);
46965d8bf0aSStefano Zampini           Stilda = sub_schurs->S_Ej[i];
47065d8bf0aSStefano Zampini         } else {
47165d8bf0aSStefano Zampini           KSP ksp;
47265d8bf0aSStefano Zampini           PC  pc;
47365d8bf0aSStefano Zampini           Mat S_EF,S_FE,S_FF,Stilda_impl;
47465d8bf0aSStefano Zampini           IS  is_F;
4752b973097SStefano Zampini           PetscScalar eps=1.e-8;
4762b973097SStefano Zampini           PetscBool chop=PETSC_FALSE;
47765d8bf0aSStefano Zampini 
47865d8bf0aSStefano Zampini           ierr = ISComplement(is_E,0,n_all,&is_F);CHKERRQ(ierr);
47965d8bf0aSStefano Zampini           ierr = MatGetSubMatrix(S_all,is_E,is_F,MAT_INITIAL_MATRIX,&S_EF);CHKERRQ(ierr);
48065d8bf0aSStefano Zampini           ierr = MatGetSubMatrix(S_all,is_F,is_F,MAT_INITIAL_MATRIX,&S_FF);CHKERRQ(ierr);
48165d8bf0aSStefano Zampini           ierr = MatGetSubMatrix(S_all,is_F,is_E,MAT_INITIAL_MATRIX,&S_FE);CHKERRQ(ierr);
482a1337663SStefano Zampini           ierr = ISDestroy(&is_F);CHKERRQ(ierr);
4832b973097SStefano Zampini           if (chop) {
4842b973097SStefano Zampini             ierr = MatChop(S_FF,eps);CHKERRQ(ierr);
4852b973097SStefano Zampini             ierr = MatConvert(S_FF,MATAIJ,MAT_REUSE_MATRIX,&S_FF);CHKERRQ(ierr);
4862b973097SStefano Zampini           }
487a1337663SStefano Zampini           ierr = MatCreateSchurComplement(S_FF,S_FF,S_FE,S_EF,sub_schurs->S_Ej[i],&Stilda_impl);CHKERRQ(ierr);
488a1337663SStefano Zampini           ierr = MatDestroy(&S_FF);CHKERRQ(ierr);
489a1337663SStefano Zampini           ierr = MatDestroy(&S_FE);CHKERRQ(ierr);
490a1337663SStefano Zampini           ierr = MatDestroy(&S_EF);CHKERRQ(ierr);
49165d8bf0aSStefano Zampini           ierr = MatSchurComplementGetKSP(Stilda_impl,&ksp);CHKERRQ(ierr);
49265d8bf0aSStefano Zampini           ierr = KSPSetType(ksp,KSPPREONLY);CHKERRQ(ierr);
49365d8bf0aSStefano Zampini           ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr);
49465d8bf0aSStefano Zampini           if (is_symmetric) {
49565d8bf0aSStefano Zampini             ierr = PCSetType(pc,PCCHOLESKY);CHKERRQ(ierr);
49665d8bf0aSStefano Zampini           } else {
49765d8bf0aSStefano Zampini             ierr = PCSetType(pc,PCLU);CHKERRQ(ierr);
49865d8bf0aSStefano Zampini           }
4992b973097SStefano Zampini           if (!chop) {
50065d8bf0aSStefano Zampini             ierr = PCFactorSetUseInPlace(pc,PETSC_TRUE);CHKERRQ(ierr);
5012b973097SStefano Zampini           } else {
5022b973097SStefano Zampini             ierr = PCFactorSetMatSolverPackage(pc,MATSOLVERMUMPS);CHKERRQ(ierr);
5032b973097SStefano Zampini           }
50465d8bf0aSStefano Zampini           ierr = KSPSetUp(ksp);CHKERRQ(ierr);
50565d8bf0aSStefano Zampini           ierr = PCBDDCComputeExplicitSchur(Stilda_impl,&Stilda);CHKERRQ(ierr);
50665d8bf0aSStefano Zampini           ierr = MatDestroy(&Stilda_impl);CHKERRQ(ierr);
50765d8bf0aSStefano Zampini         }
5082b973097SStefano Zampini /*
5092b973097SStefano Zampini         PetscViewerPushFormat(PETSC_VIEWER_STDOUT_SELF,PETSC_VIEWER_ASCII_MATLAB);
5102b973097SStefano Zampini         ierr = MatView(Stilda,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
5112b973097SStefano Zampini         PetscViewerPopFormat(PETSC_VIEWER_STDOUT_SELF);
5122b973097SStefano Zampini */
513b96c3477SStefano Zampini         ierr = MatDestroy(&sub_schurs->S_Ej_tilda[i]);CHKERRQ(ierr);
5142b973097SStefano Zampini         sub_schurs->S_Ej_tilda[i] = Stilda;
51508122e43SStefano Zampini 
51608122e43SStefano Zampini 
517a1337663SStefano Zampini         ierr = MatDenseGetArray(sub_schurs->S_Ej_tilda[i],&vals);CHKERRQ(ierr);
51808122e43SStefano Zampini         if (deluxe) { /* when using deluxe scaling, we need (S_1^-1+S_2^-1)^-1 */
51908122e43SStefano Zampini           PetscScalar *vals2;
52008122e43SStefano Zampini 
52108122e43SStefano Zampini           ierr = MatDenseGetArray(sub_schurs->S_Ej[i],&vals2);CHKERRQ(ierr);
52208122e43SStefano Zampini           /* need to be optimized (cholesky) */
52308122e43SStefano Zampini           ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr);
52408122e43SStefano Zampini           ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
52508122e43SStefano Zampini           PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,vals,&B_N,pivots,&B_ierr));
52608122e43SStefano Zampini           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
52708122e43SStefano Zampini           PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,vals,&B_N,pivots,work,&B_lwork,&B_ierr));
52808122e43SStefano Zampini           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
52908122e43SStefano Zampini           PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,vals2,&B_N,pivots,&B_ierr));
53008122e43SStefano Zampini           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
53108122e43SStefano Zampini           PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,vals2,&B_N,pivots,work,&B_lwork,&B_ierr));
53208122e43SStefano Zampini           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
53308122e43SStefano Zampini           ierr = PetscFPTrapPop();CHKERRQ(ierr);
53408122e43SStefano Zampini           ierr = MatSetValues(S_Ej_inv_all,subset_size,dummy_idx,subset_size,dummy_idx,vals2,INSERT_VALUES);CHKERRQ(ierr);
53508122e43SStefano Zampini           ierr = MatDenseRestoreArray(sub_schurs->S_Ej[i],&vals2);CHKERRQ(ierr);
53608122e43SStefano Zampini         }
537a1337663SStefano Zampini         ierr = MatSetValues(S_Ej_tilda_all,subset_size,dummy_idx,subset_size,dummy_idx,vals,INSERT_VALUES);CHKERRQ(ierr);
538a1337663SStefano Zampini         ierr = MatDenseRestoreArray(sub_schurs->S_Ej_tilda[i],&vals);CHKERRQ(ierr);
5392b973097SStefano Zampini       } else {
540b96c3477SStefano Zampini         ierr = MatDestroy(&sub_schurs->S_Ej_tilda[i]);CHKERRQ(ierr);
54165d8bf0aSStefano Zampini       }
54265d8bf0aSStefano Zampini       ierr = ISDestroy(&is_E);CHKERRQ(ierr);
543883469d8SStefano Zampini       local_size += subset_size;
544883469d8SStefano Zampini     }
545883469d8SStefano Zampini     ierr = PetscFree(dummy_idx);CHKERRQ(ierr);
5465db18549SStefano Zampini   }
547a1337663SStefano Zampini   ierr = PetscFree(nnz);CHKERRQ(ierr);
548a1337663SStefano Zampini   ierr = MatDestroy(&S_all);CHKERRQ(ierr);
5495db18549SStefano Zampini   ierr = MatAssemblyBegin(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
5505db18549SStefano Zampini   ierr = MatAssemblyEnd(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
551a1337663SStefano Zampini   if (compute_Stilda) {
552a1337663SStefano Zampini     ierr = MatAssemblyBegin(S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
553a1337663SStefano Zampini     ierr = MatAssemblyEnd(S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
55408122e43SStefano Zampini     if (deluxe) {
55508122e43SStefano Zampini       ierr = MatAssemblyBegin(S_Ej_inv_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
55608122e43SStefano Zampini       ierr = MatAssemblyEnd(S_Ej_inv_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
55708122e43SStefano Zampini     }
558a1337663SStefano Zampini   }
559a1337663SStefano Zampini 
5605db18549SStefano Zampini   /* Global matrix of all assembled Schur on subsets */
561883469d8SStefano 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);
5625db18549SStefano Zampini   ierr = ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)sub_schurs->l2gmap),1,local_size,all_local_idx_G,PETSC_COPY_VALUES,&l2gmap_subsets);CHKERRQ(ierr);
5635db18549SStefano Zampini   ierr = MatCreateIS(PetscObjectComm((PetscObject)sub_schurs->l2gmap),1,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size,l2gmap_subsets,&work_mat);CHKERRQ(ierr);
5645db18549SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&l2gmap_subsets);CHKERRQ(ierr);
5655db18549SStefano Zampini   ierr = MatISSetLocalMat(work_mat,sub_schurs->S_Ej_all);CHKERRQ(ierr);
566*3927de2eSStefano Zampini   ierr = MatCreate(PetscObjectComm((PetscObject)work_mat),&global_schur_subsets);CHKERRQ(ierr);
567*3927de2eSStefano Zampini   ierr = MatSetSizes(global_schur_subsets,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size);CHKERRQ(ierr);
568*3927de2eSStefano Zampini   ierr = MatSetType(global_schur_subsets,MATMPIAIJ);CHKERRQ(ierr);
569*3927de2eSStefano Zampini   ierr = MatISSetMPIXAIJPreallocation_Private(work_mat,global_schur_subsets,PETSC_TRUE);CHKERRQ(ierr);
570*3927de2eSStefano Zampini   ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
5715db18549SStefano Zampini   /* Get local part of (\sum_j S_Ej) */
572883469d8SStefano Zampini   ierr = PetscFree(all_local_idx_N);CHKERRQ(ierr);
573883469d8SStefano Zampini   ierr = ISCreateGeneral(PetscObjectComm((PetscObject)sub_schurs->l2gmap),local_size,all_local_idx_G,PETSC_OWN_POINTER,&temp_is);CHKERRQ(ierr);
574b96c3477SStefano Zampini   if (sub_schurs->sum_S_Ej_all) {
575b96c3477SStefano Zampini     ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,MAT_REUSE_MATRIX,&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
576b96c3477SStefano Zampini   } else {
577b96c3477SStefano Zampini     ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,MAT_INITIAL_MATRIX,&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
578b96c3477SStefano Zampini   }
57908122e43SStefano Zampini 
58008122e43SStefano Zampini   if (Stilda_computed) {
581a1337663SStefano Zampini     ierr = MatISSetLocalMat(work_mat,S_Ej_tilda_all);CHKERRQ(ierr);
582a1337663SStefano Zampini     ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
583b96c3477SStefano Zampini     if (sub_schurs->sum_S_Ej_tilda_all) {
584b96c3477SStefano Zampini       ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,MAT_REUSE_MATRIX,&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
585b96c3477SStefano Zampini     } else {
586b96c3477SStefano Zampini       ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,MAT_INITIAL_MATRIX,&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
587b96c3477SStefano Zampini     }
58808122e43SStefano Zampini     if (deluxe) {
58908122e43SStefano Zampini       PetscInt    cum;
59008122e43SStefano Zampini       PetscScalar *array,*array2;
59108122e43SStefano Zampini       ierr = MatISSetLocalMat(work_mat,S_Ej_inv_all);CHKERRQ(ierr);
59208122e43SStefano Zampini       ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr);
59308122e43SStefano Zampini       if (sub_schurs->sum_S_Ej_inv_all) {
59408122e43SStefano Zampini         ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,MAT_REUSE_MATRIX,&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
59508122e43SStefano Zampini       } else {
59608122e43SStefano Zampini         ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,MAT_INITIAL_MATRIX,&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
59708122e43SStefano Zampini       }
59808122e43SStefano Zampini       /* invert blocks of sum_S_Ej_inv_all */
59908122e43SStefano Zampini       ierr = MatSeqAIJGetArray(sub_schurs->sum_S_Ej_inv_all,&array);CHKERRQ(ierr);
60008122e43SStefano Zampini       ierr = MatSeqAIJGetArray(sub_schurs->sum_S_Ej_tilda_all,&array2);CHKERRQ(ierr);
60108122e43SStefano Zampini       cum = 0;
60208122e43SStefano Zampini       for (i=0;i<sub_schurs->n_subs;i++) {
60308122e43SStefano Zampini         ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr);
60408122e43SStefano Zampini         if (PetscBTLookup(sub_schurs->computed_Stilda_subs,i)) {
60508122e43SStefano Zampini           /* need to be optimized (cholesky) */
60608122e43SStefano Zampini           ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr);
60708122e43SStefano Zampini           ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
60808122e43SStefano Zampini           PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,array+cum,&B_N,pivots,&B_ierr));
60908122e43SStefano Zampini           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
61008122e43SStefano Zampini           PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,array+cum,&B_N,pivots,work,&B_lwork,&B_ierr));
61108122e43SStefano Zampini           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
61208122e43SStefano Zampini           PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,array2+cum,&B_N,pivots,&B_ierr));
61308122e43SStefano Zampini           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
61408122e43SStefano Zampini           PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,array2+cum,&B_N,pivots,work,&B_lwork,&B_ierr));
61508122e43SStefano Zampini           if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
61608122e43SStefano Zampini           ierr = PetscFPTrapPop();CHKERRQ(ierr);
61708122e43SStefano Zampini         }
61808122e43SStefano Zampini         cum += subset_size*subset_size;
61908122e43SStefano Zampini       }
62008122e43SStefano Zampini       ierr = MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_inv_all,&array);CHKERRQ(ierr);
62108122e43SStefano Zampini       ierr = MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_tilda_all,&array2);CHKERRQ(ierr);
62208122e43SStefano Zampini     }
6233202ece2SStefano Zampini   }
6243202ece2SStefano Zampini 
62508122e43SStefano Zampini   ierr = PetscFree2(work,pivots);CHKERRQ(ierr);
626a1337663SStefano Zampini   ierr = MatDestroy(&global_schur_subsets);CHKERRQ(ierr);
627a1337663SStefano Zampini   ierr = MatDestroy(&S_Ej_tilda_all);CHKERRQ(ierr);
62808122e43SStefano Zampini   ierr = MatDestroy(&S_Ej_inv_all);CHKERRQ(ierr);
6293202ece2SStefano Zampini   ierr = MatDestroy(&work_mat);CHKERRQ(ierr);
6305db18549SStefano Zampini   ierr = ISDestroy(&temp_is);CHKERRQ(ierr);
631b1b3d7a2SStefano Zampini   PetscFunctionReturn(0);
632b1b3d7a2SStefano Zampini }
633b1b3d7a2SStefano Zampini 
634b1b3d7a2SStefano Zampini #undef __FUNCT__
635b1b3d7a2SStefano Zampini #define __FUNCT__ "PCBDDCSubSchursInit"
6365db18549SStefano Zampini PetscErrorCode PCBDDCSubSchursInit(PCBDDCSubSchurs sub_schurs, Mat A, Mat S, IS is_I, IS is_B, PCBDDCGraph graph, ISLocalToGlobalMapping BtoNmap, PetscInt seqthreshold)
637b1b3d7a2SStefano Zampini {
6389bb4a8caSStefano Zampini   IS                  *faces,*edges,*all_cc,vertices;
639b1b3d7a2SStefano Zampini   PetscInt            *index_sequential,*index_parallel;
640b1b3d7a2SStefano Zampini   PetscInt            *auxlocal_sequential,*auxlocal_parallel;
641b1b3d7a2SStefano Zampini   PetscInt            *auxglobal_sequential,*auxglobal_parallel;
642b96c3477SStefano Zampini   PetscInt            *auxmapping;
643b1b3d7a2SStefano Zampini   PetscInt            i,max_subset_size;
644b1b3d7a2SStefano Zampini   PetscInt            n_sequential_problems,n_local_sequential_problems,n_parallel_problems,n_local_parallel_problems;
645b1b3d7a2SStefano Zampini   PetscInt            n_faces,n_edges,n_all_cc;
646b1b3d7a2SStefano Zampini   PetscBool           is_sorted;
647b1b3d7a2SStefano Zampini   PetscErrorCode  ierr;
648b1b3d7a2SStefano Zampini 
649b1b3d7a2SStefano Zampini   PetscFunctionBegin;
650b1b3d7a2SStefano Zampini   ierr = ISSorted(is_I,&is_sorted);CHKERRQ(ierr);
651b1b3d7a2SStefano Zampini   if (!is_sorted) {
652b1b3d7a2SStefano Zampini     SETERRQ(PetscObjectComm((PetscObject)is_I),PETSC_ERR_PLIB,"IS for I dofs should be shorted");
653b1b3d7a2SStefano Zampini   }
654b1b3d7a2SStefano Zampini   ierr = ISSorted(is_B,&is_sorted);CHKERRQ(ierr);
655b1b3d7a2SStefano Zampini   if (!is_sorted) {
656b1b3d7a2SStefano Zampini     SETERRQ(PetscObjectComm((PetscObject)is_B),PETSC_ERR_PLIB,"IS for B dofs should be shorted");
657b1b3d7a2SStefano Zampini   }
658b1b3d7a2SStefano Zampini 
659b1b3d7a2SStefano Zampini   /* reset any previous data */
660b1b3d7a2SStefano Zampini   ierr = PCBDDCSubSchursReset(sub_schurs);CHKERRQ(ierr);
661b1b3d7a2SStefano Zampini 
662b1b3d7a2SStefano Zampini   /* get index sets for faces and edges */
6639bb4a8caSStefano Zampini   ierr = PCBDDCGraphGetCandidatesIS(graph,&n_faces,&faces,&n_edges,&edges,&vertices);CHKERRQ(ierr);
664b1b3d7a2SStefano Zampini   n_all_cc = n_faces+n_edges;
66508122e43SStefano Zampini   ierr = PetscBTCreate(n_all_cc,&sub_schurs->is_edge);CHKERRQ(ierr);
66608122e43SStefano Zampini   ierr = PetscBTCreate(n_all_cc,&sub_schurs->computed_Stilda_subs);CHKERRQ(ierr);
667b1b3d7a2SStefano Zampini   ierr = PetscMalloc1(n_all_cc,&all_cc);CHKERRQ(ierr);
668b1b3d7a2SStefano Zampini   for (i=0;i<n_faces;i++) {
669b1b3d7a2SStefano Zampini     all_cc[i] = faces[i];
670b1b3d7a2SStefano Zampini   }
671b1b3d7a2SStefano Zampini   for (i=0;i<n_edges;i++) {
672b1b3d7a2SStefano Zampini     all_cc[n_faces+i] = edges[i];
67308122e43SStefano Zampini     ierr = PetscBTSet(sub_schurs->is_edge,n_faces+i);CHKERRQ(ierr);
674b1b3d7a2SStefano Zampini   }
675b1b3d7a2SStefano Zampini   ierr = PetscFree(faces);CHKERRQ(ierr);
676b1b3d7a2SStefano Zampini   ierr = PetscFree(edges);CHKERRQ(ierr);
677b1b3d7a2SStefano Zampini 
678b1b3d7a2SStefano Zampini   /* map interface's subsets */
679b1b3d7a2SStefano Zampini   max_subset_size = 0;
680b1b3d7a2SStefano Zampini   for (i=0;i<n_all_cc;i++) {
681b1b3d7a2SStefano Zampini     PetscInt subset_size;
682b1b3d7a2SStefano Zampini     ierr = ISGetLocalSize(all_cc[i],&subset_size);CHKERRQ(ierr);
683b1b3d7a2SStefano Zampini     max_subset_size = PetscMax(max_subset_size,subset_size);
684b1b3d7a2SStefano Zampini   }
685b1b3d7a2SStefano Zampini   ierr = PetscMalloc1(max_subset_size,&auxmapping);CHKERRQ(ierr);
686b1b3d7a2SStefano Zampini   ierr = PetscMalloc2(graph->ncc,&auxlocal_sequential,
687b1b3d7a2SStefano Zampini                       graph->ncc,&auxlocal_parallel);CHKERRQ(ierr);
688b1b3d7a2SStefano Zampini   ierr = PetscMalloc2(graph->ncc,&index_sequential,
689b1b3d7a2SStefano Zampini                       graph->ncc,&index_parallel);CHKERRQ(ierr);
690b1b3d7a2SStefano Zampini 
691883469d8SStefano Zampini   /* if threshold is negative uses all sequential problems (possibly using MUMPS) */
692883469d8SStefano Zampini   sub_schurs->use_mumps = PETSC_FALSE;
693883469d8SStefano Zampini   if (seqthreshold < 0) {
694883469d8SStefano Zampini     seqthreshold = max_subset_size;
695883469d8SStefano Zampini #if defined(PETSC_HAVE_MUMPS)
696883469d8SStefano Zampini     sub_schurs->use_mumps = !!A;
697883469d8SStefano Zampini #endif
698883469d8SStefano Zampini   }
699883469d8SStefano Zampini 
700b1b3d7a2SStefano Zampini 
701b1b3d7a2SStefano Zampini   /* determine which problem has to be solved in parallel or sequentially */
702b1b3d7a2SStefano Zampini   n_local_sequential_problems = 0;
703b1b3d7a2SStefano Zampini   n_local_parallel_problems = 0;
704b1b3d7a2SStefano Zampini   for (i=0;i<n_all_cc;i++) {
705b1b3d7a2SStefano Zampini     PetscInt       subset_size,j,min_loc = 0;
706b1b3d7a2SStefano Zampini     const PetscInt *idxs;
707b1b3d7a2SStefano Zampini 
708b1b3d7a2SStefano Zampini     ierr = ISGetLocalSize(all_cc[i],&subset_size);CHKERRQ(ierr);
709b1b3d7a2SStefano Zampini     ierr = ISGetIndices(all_cc[i],&idxs);CHKERRQ(ierr);
710b1b3d7a2SStefano Zampini     ierr = ISLocalToGlobalMappingApply(graph->l2gmap,subset_size,idxs,auxmapping);CHKERRQ(ierr);
711b1b3d7a2SStefano Zampini     for (j=1;j<subset_size;j++) {
712b1b3d7a2SStefano Zampini       if (auxmapping[j]<auxmapping[min_loc]) {
713b1b3d7a2SStefano Zampini         min_loc = j;
714b1b3d7a2SStefano Zampini       }
715b1b3d7a2SStefano Zampini     }
716b1b3d7a2SStefano Zampini     if (subset_size > seqthreshold) {
717b1b3d7a2SStefano Zampini       index_parallel[n_local_parallel_problems] = i;
718b1b3d7a2SStefano Zampini       auxlocal_parallel[n_local_parallel_problems] = idxs[min_loc];
719b1b3d7a2SStefano Zampini       n_local_parallel_problems++;
720b1b3d7a2SStefano Zampini     } else {
721b1b3d7a2SStefano Zampini       index_sequential[n_local_sequential_problems] = i;
722b1b3d7a2SStefano Zampini       auxlocal_sequential[n_local_sequential_problems] = idxs[min_loc];
723b1b3d7a2SStefano Zampini       n_local_sequential_problems++;
724b1b3d7a2SStefano Zampini     }
725b1b3d7a2SStefano Zampini     ierr = ISRestoreIndices(all_cc[i],&idxs);CHKERRQ(ierr);
726b1b3d7a2SStefano Zampini   }
727b1b3d7a2SStefano Zampini 
728b1b3d7a2SStefano Zampini   /* Number parallel problems */
729b1b3d7a2SStefano Zampini   auxglobal_parallel = 0;
730b1b3d7a2SStefano Zampini   ierr = PCBDDCSubsetNumbering(PetscObjectComm((PetscObject)graph->l2gmap),graph->l2gmap,n_local_parallel_problems,auxlocal_parallel,PETSC_NULL,&n_parallel_problems,&auxglobal_parallel);CHKERRQ(ierr);
731b1b3d7a2SStefano Zampini 
732b1b3d7a2SStefano Zampini   /* Number sequential problems */
733b1b3d7a2SStefano Zampini   auxglobal_sequential = 0;
734b1b3d7a2SStefano Zampini   ierr = PCBDDCSubsetNumbering(PetscObjectComm((PetscObject)graph->l2gmap),graph->l2gmap,n_local_sequential_problems,auxlocal_sequential,PETSC_NULL,&n_sequential_problems,&auxglobal_sequential);CHKERRQ(ierr);
735b1b3d7a2SStefano Zampini 
736b1b3d7a2SStefano Zampini   /* update info in sub_schurs */
737883469d8SStefano Zampini   if (sub_schurs->use_mumps && A) {
7381e9c79c2SStefano Zampini     ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
7391e9c79c2SStefano Zampini     sub_schurs->A = A;
7401e9c79c2SStefano Zampini   }
741b1b3d7a2SStefano Zampini   ierr = PetscObjectReference((PetscObject)S);CHKERRQ(ierr);
742b1b3d7a2SStefano Zampini   sub_schurs->S = S;
743b1b3d7a2SStefano Zampini   ierr = PetscObjectReference((PetscObject)is_I);CHKERRQ(ierr);
744b1b3d7a2SStefano Zampini   sub_schurs->is_I = is_I;
745b1b3d7a2SStefano Zampini   ierr = PetscObjectReference((PetscObject)is_B);CHKERRQ(ierr);
746b1b3d7a2SStefano Zampini   sub_schurs->is_B = is_B;
7475db18549SStefano Zampini   ierr = PetscObjectReference((PetscObject)graph->l2gmap);CHKERRQ(ierr);
7485db18549SStefano Zampini   sub_schurs->l2gmap = graph->l2gmap;
7495db18549SStefano Zampini   ierr = PetscObjectReference((PetscObject)BtoNmap);CHKERRQ(ierr);
7505db18549SStefano Zampini   sub_schurs->BtoNmap = BtoNmap;
751b1b3d7a2SStefano Zampini   sub_schurs->n_subs_seq = n_local_sequential_problems;
752b1b3d7a2SStefano Zampini   sub_schurs->n_subs_par = n_local_parallel_problems;
753b1b3d7a2SStefano Zampini   sub_schurs->n_subs_seq_g = n_sequential_problems;
754b1b3d7a2SStefano Zampini   sub_schurs->n_subs_par_g = n_parallel_problems;
755b1b3d7a2SStefano Zampini   sub_schurs->n_subs = sub_schurs->n_subs_seq + sub_schurs->n_subs_par;
756b1b3d7a2SStefano Zampini   sub_schurs->is_subs = all_cc;
7579bb4a8caSStefano Zampini   if (!sub_schurs->use_mumps) { /* for adaptive selection */
758b96c3477SStefano Zampini     for (i=0;i<sub_schurs->n_subs;i++) {
759b96c3477SStefano Zampini       ierr = ISSort(sub_schurs->is_subs[i]);CHKERRQ(ierr);
760b96c3477SStefano Zampini     }
7619bb4a8caSStefano Zampini   }
7629bb4a8caSStefano Zampini   sub_schurs->is_Ej_com = vertices;
763b1b3d7a2SStefano Zampini   sub_schurs->index_sequential = index_sequential;
764b1b3d7a2SStefano Zampini   sub_schurs->index_parallel = index_parallel;
765b1b3d7a2SStefano Zampini   sub_schurs->auxglobal_sequential = auxglobal_sequential;
766b1b3d7a2SStefano Zampini   sub_schurs->auxglobal_parallel = auxglobal_parallel;
767b1b3d7a2SStefano Zampini 
768b96c3477SStefano Zampini 
769b96c3477SStefano Zampini   /* allocate space for schur complements */
770b96c3477SStefano Zampini   ierr = PetscCalloc2(sub_schurs->n_subs,&sub_schurs->S_Ej,
771b96c3477SStefano Zampini                       sub_schurs->n_subs,&sub_schurs->S_Ej_tilda);CHKERRQ(ierr);
772b96c3477SStefano Zampini   sub_schurs->S_Ej_all = NULL;
773b96c3477SStefano Zampini   sub_schurs->sum_S_Ej_all = NULL;
77408122e43SStefano Zampini   sub_schurs->sum_S_Ej_inv_all = NULL;
775b96c3477SStefano Zampini   sub_schurs->sum_S_Ej_tilda_all = NULL;
776b96c3477SStefano Zampini   sub_schurs->is_Ej_all = NULL;
777b96c3477SStefano Zampini 
778b1b3d7a2SStefano Zampini   /* free workspace */
779b1b3d7a2SStefano Zampini   ierr = PetscFree(auxmapping);CHKERRQ(ierr);
780b1b3d7a2SStefano Zampini   ierr = PetscFree2(auxlocal_sequential,auxlocal_parallel);CHKERRQ(ierr);
781b1b3d7a2SStefano Zampini 
782b1b3d7a2SStefano Zampini   PetscFunctionReturn(0);
783b1b3d7a2SStefano Zampini }
784b1b3d7a2SStefano Zampini 
785b1b3d7a2SStefano Zampini #undef __FUNCT__
78634a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCSubSchursCreate"
78734a97f8cSStefano Zampini PetscErrorCode PCBDDCSubSchursCreate(PCBDDCSubSchurs *sub_schurs)
78834a97f8cSStefano Zampini {
78934a97f8cSStefano Zampini   PCBDDCSubSchurs schurs_ctx;
79034a97f8cSStefano Zampini   PetscErrorCode  ierr;
79134a97f8cSStefano Zampini 
79234a97f8cSStefano Zampini   PetscFunctionBegin;
79334a97f8cSStefano Zampini   ierr = PetscNew(&schurs_ctx);CHKERRQ(ierr);
7945ff63025SStefano Zampini   schurs_ctx->n_subs = 0;
79534a97f8cSStefano Zampini   *sub_schurs = schurs_ctx;
79634a97f8cSStefano Zampini   PetscFunctionReturn(0);
79734a97f8cSStefano Zampini }
79834a97f8cSStefano Zampini 
79934a97f8cSStefano Zampini #undef __FUNCT__
80034a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCSubSchursDestroy"
80134a97f8cSStefano Zampini PetscErrorCode PCBDDCSubSchursDestroy(PCBDDCSubSchurs *sub_schurs)
80234a97f8cSStefano Zampini {
80334a97f8cSStefano Zampini   PetscErrorCode ierr;
80434a97f8cSStefano Zampini 
80534a97f8cSStefano Zampini   PetscFunctionBegin;
80634a97f8cSStefano Zampini   ierr = PCBDDCSubSchursReset(*sub_schurs);CHKERRQ(ierr);
80734a97f8cSStefano Zampini   ierr = PetscFree(*sub_schurs);CHKERRQ(ierr);
80834a97f8cSStefano Zampini   PetscFunctionReturn(0);
80934a97f8cSStefano Zampini }
81034a97f8cSStefano Zampini 
81134a97f8cSStefano Zampini #undef __FUNCT__
81234a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCSubSchursReset"
81334a97f8cSStefano Zampini PetscErrorCode PCBDDCSubSchursReset(PCBDDCSubSchurs sub_schurs)
81434a97f8cSStefano Zampini {
81534a97f8cSStefano Zampini   PetscInt       i;
81634a97f8cSStefano Zampini   PetscErrorCode ierr;
81734a97f8cSStefano Zampini 
81834a97f8cSStefano Zampini   PetscFunctionBegin;
8191e9c79c2SStefano Zampini   ierr = MatDestroy(&sub_schurs->A);CHKERRQ(ierr);
820b1b3d7a2SStefano Zampini   ierr = MatDestroy(&sub_schurs->S);CHKERRQ(ierr);
821b1b3d7a2SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_I);CHKERRQ(ierr);
822b1b3d7a2SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_B);CHKERRQ(ierr);
8235db18549SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->l2gmap);CHKERRQ(ierr);
8245db18549SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->BtoNmap);CHKERRQ(ierr);
82541c3ba1bSStefano Zampini   ierr = MatDestroy(&sub_schurs->S_Ej_all);CHKERRQ(ierr);
82641c3ba1bSStefano Zampini   ierr = MatDestroy(&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr);
82708122e43SStefano Zampini   ierr = MatDestroy(&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr);
828a1337663SStefano Zampini   ierr = MatDestroy(&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr);
8295db18549SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_Ej_all);CHKERRQ(ierr);
8309bb4a8caSStefano Zampini   ierr = ISDestroy(&sub_schurs->is_Ej_com);CHKERRQ(ierr);
83108122e43SStefano Zampini   ierr = PetscBTDestroy(&sub_schurs->is_edge);CHKERRQ(ierr);
83208122e43SStefano Zampini   ierr = PetscBTDestroy(&sub_schurs->computed_Stilda_subs);CHKERRQ(ierr);
83334a97f8cSStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
834b1b3d7a2SStefano Zampini     ierr = ISDestroy(&sub_schurs->is_subs[i]);CHKERRQ(ierr);
83534a97f8cSStefano Zampini     ierr = MatDestroy(&sub_schurs->S_Ej[i]);CHKERRQ(ierr);
8362b973097SStefano Zampini     ierr = MatDestroy(&sub_schurs->S_Ej_tilda[i]);CHKERRQ(ierr);
83734a97f8cSStefano Zampini   }
83868270318SStefano Zampini   ierr = ISDestroy(&sub_schurs->is_I_layer);CHKERRQ(ierr);
8395ff63025SStefano Zampini   if (sub_schurs->n_subs) {
840b1b3d7a2SStefano Zampini     ierr = PetscFree(sub_schurs->is_subs);CHKERRQ(ierr);
841b96c3477SStefano Zampini     ierr = PetscFree2(sub_schurs->S_Ej,sub_schurs->S_Ej_tilda);CHKERRQ(ierr);
842b1b3d7a2SStefano Zampini     ierr = PetscFree2(sub_schurs->index_sequential,sub_schurs->index_parallel);CHKERRQ(ierr);
843b1b3d7a2SStefano Zampini     ierr = PetscFree(sub_schurs->auxglobal_sequential);CHKERRQ(ierr);
844b1b3d7a2SStefano Zampini     ierr = PetscFree(sub_schurs->auxglobal_parallel);CHKERRQ(ierr);
8455ff63025SStefano Zampini   }
84634a97f8cSStefano Zampini   sub_schurs->n_subs = 0;
84734a97f8cSStefano Zampini   PetscFunctionReturn(0);
84834a97f8cSStefano Zampini }
84934a97f8cSStefano Zampini 
85034a97f8cSStefano Zampini #undef __FUNCT__
85134a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCAdjGetNextLayer_Private"
8522a155e38SStefano Zampini PETSC_STATIC_INLINE PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt* queue_tip,PetscInt n_prev,PetscBT touched,PetscInt* xadj,PetscInt* adjncy,PetscInt* n_added)
85334a97f8cSStefano Zampini {
85434a97f8cSStefano Zampini   PetscInt       i,j,n;
85534a97f8cSStefano Zampini   PetscErrorCode ierr;
85634a97f8cSStefano Zampini 
85734a97f8cSStefano Zampini   PetscFunctionBegin;
85834a97f8cSStefano Zampini   n = 0;
85934a97f8cSStefano Zampini   for (i=-n_prev;i<0;i++) {
86034a97f8cSStefano Zampini     PetscInt start_dof = queue_tip[i];
86134a97f8cSStefano Zampini     for (j=xadj[start_dof];j<xadj[start_dof+1];j++) {
86234a97f8cSStefano Zampini       PetscInt dof = adjncy[j];
86334a97f8cSStefano Zampini       if (!PetscBTLookup(touched,dof)) {
86434a97f8cSStefano Zampini         ierr = PetscBTSet(touched,dof);CHKERRQ(ierr);
86534a97f8cSStefano Zampini         queue_tip[n] = dof;
86634a97f8cSStefano Zampini         n++;
86734a97f8cSStefano Zampini       }
86834a97f8cSStefano Zampini     }
86934a97f8cSStefano Zampini   }
87034a97f8cSStefano Zampini   *n_added = n;
87134a97f8cSStefano Zampini   PetscFunctionReturn(0);
87234a97f8cSStefano Zampini }
873