xref: /petsc/src/ksp/pc/impls/bddc/bddcschurs.c (revision 34a97f8c781d28c390054b7fe2045ae8e3fc2469)
1*34a97f8cSStefano Zampini #include <../src/ksp/pc/impls/bddc/bddc.h>
2*34a97f8cSStefano Zampini #include <../src/ksp/pc/impls/bddc/bddcprivate.h>
3*34a97f8cSStefano Zampini 
4*34a97f8cSStefano Zampini static PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt*,PetscInt,PetscBT,PetscInt*,PetscInt*,PetscInt*);
5*34a97f8cSStefano Zampini 
6*34a97f8cSStefano Zampini #undef __FUNCT__
7*34a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCSubSchursCreate"
8*34a97f8cSStefano Zampini PetscErrorCode PCBDDCSubSchursCreate(PCBDDCSubSchurs *sub_schurs)
9*34a97f8cSStefano Zampini {
10*34a97f8cSStefano Zampini   PCBDDCSubSchurs schurs_ctx;
11*34a97f8cSStefano Zampini   PetscErrorCode  ierr;
12*34a97f8cSStefano Zampini 
13*34a97f8cSStefano Zampini   PetscFunctionBegin;
14*34a97f8cSStefano Zampini   ierr = PetscNew(&schurs_ctx);CHKERRQ(ierr);
15*34a97f8cSStefano Zampini   *sub_schurs = schurs_ctx;
16*34a97f8cSStefano Zampini   PetscFunctionReturn(0);
17*34a97f8cSStefano Zampini }
18*34a97f8cSStefano Zampini 
19*34a97f8cSStefano Zampini #undef __FUNCT__
20*34a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCSubSchursDestroy"
21*34a97f8cSStefano Zampini PetscErrorCode PCBDDCSubSchursDestroy(PCBDDCSubSchurs *sub_schurs)
22*34a97f8cSStefano Zampini {
23*34a97f8cSStefano Zampini   PetscErrorCode ierr;
24*34a97f8cSStefano Zampini 
25*34a97f8cSStefano Zampini   PetscFunctionBegin;
26*34a97f8cSStefano Zampini   ierr = PCBDDCSubSchursReset(*sub_schurs);CHKERRQ(ierr);
27*34a97f8cSStefano Zampini   ierr = PetscFree(*sub_schurs);CHKERRQ(ierr);
28*34a97f8cSStefano Zampini   PetscFunctionReturn(0);
29*34a97f8cSStefano Zampini }
30*34a97f8cSStefano Zampini 
31*34a97f8cSStefano Zampini #undef __FUNCT__
32*34a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCSubSchursReset"
33*34a97f8cSStefano Zampini PetscErrorCode PCBDDCSubSchursReset(PCBDDCSubSchurs sub_schurs)
34*34a97f8cSStefano Zampini {
35*34a97f8cSStefano Zampini   PetscInt       i;
36*34a97f8cSStefano Zampini   PetscErrorCode ierr;
37*34a97f8cSStefano Zampini 
38*34a97f8cSStefano Zampini   PetscFunctionBegin;
39*34a97f8cSStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
40*34a97f8cSStefano Zampini     ierr = ISDestroy(&sub_schurs->is_AEj_I[i]);CHKERRQ(ierr);
41*34a97f8cSStefano Zampini     ierr = ISDestroy(&sub_schurs->is_AEj_B[i]);CHKERRQ(ierr);
42*34a97f8cSStefano Zampini     ierr = MatDestroy(&sub_schurs->S_Ej[i]);CHKERRQ(ierr);
43*34a97f8cSStefano Zampini     ierr = VecDestroy(&sub_schurs->work1[i]);CHKERRQ(ierr);
44*34a97f8cSStefano Zampini     ierr = VecDestroy(&sub_schurs->work2[i]);CHKERRQ(ierr);
45*34a97f8cSStefano Zampini   }
46*34a97f8cSStefano Zampini   ierr = PetscFree5(sub_schurs->is_AEj_I,sub_schurs->is_AEj_B,sub_schurs->S_Ej,sub_schurs->work1,sub_schurs->work2);CHKERRQ(ierr);
47*34a97f8cSStefano Zampini   sub_schurs->n_subs = 0;
48*34a97f8cSStefano Zampini   PetscFunctionReturn(0);
49*34a97f8cSStefano Zampini }
50*34a97f8cSStefano Zampini 
51*34a97f8cSStefano Zampini #undef __FUNCT__
52*34a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCAdjGetNextLayer_Private"
53*34a97f8cSStefano Zampini static PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt* queue_tip,PetscInt n_prev,PetscBT touched,PetscInt* xadj,PetscInt* adjncy,PetscInt* n_added)
54*34a97f8cSStefano Zampini {
55*34a97f8cSStefano Zampini   PetscInt       i,j,n;
56*34a97f8cSStefano Zampini   PetscErrorCode ierr;
57*34a97f8cSStefano Zampini 
58*34a97f8cSStefano Zampini   PetscFunctionBegin;
59*34a97f8cSStefano Zampini   n = 0;
60*34a97f8cSStefano Zampini   for (i=-n_prev;i<0;i++) {
61*34a97f8cSStefano Zampini     PetscInt start_dof = queue_tip[i];
62*34a97f8cSStefano Zampini     for (j=xadj[start_dof];j<xadj[start_dof+1];j++) {
63*34a97f8cSStefano Zampini       PetscInt dof = adjncy[j];
64*34a97f8cSStefano Zampini       if (!PetscBTLookup(touched,dof)) {
65*34a97f8cSStefano Zampini         ierr = PetscBTSet(touched,dof);CHKERRQ(ierr);
66*34a97f8cSStefano Zampini         queue_tip[n] = dof;
67*34a97f8cSStefano Zampini         n++;
68*34a97f8cSStefano Zampini       }
69*34a97f8cSStefano Zampini     }
70*34a97f8cSStefano Zampini   }
71*34a97f8cSStefano Zampini   *n_added = n;
72*34a97f8cSStefano Zampini   PetscFunctionReturn(0);
73*34a97f8cSStefano Zampini }
74*34a97f8cSStefano Zampini 
75*34a97f8cSStefano Zampini #undef __FUNCT__
76*34a97f8cSStefano Zampini #define __FUNCT__ "PCBDDCSubSchursSetUp"
77*34a97f8cSStefano Zampini PetscErrorCode PCBDDCSubSchursSetUp(PCBDDCSubSchurs sub_schurs, Mat S, IS is_A_I, IS is_A_B, PetscInt ncc, PetscInt cptr[], PetscInt queue[], PetscInt xadj[], PetscInt adjncy[], PetscBool skip_corners, PetscInt nlayers)
78*34a97f8cSStefano Zampini {
79*34a97f8cSStefano Zampini   Mat                    A_II,A_IB,A_BI,A_BB;
80*34a97f8cSStefano Zampini   ISLocalToGlobalMapping BtoNmap,ItoNmap;
81*34a97f8cSStefano Zampini   PetscBT                touched;
82*34a97f8cSStefano Zampini   PetscInt               i,n_I,n_B,n_local,*local_numbering;
83*34a97f8cSStefano Zampini   PetscBool              is_sorted;
84*34a97f8cSStefano Zampini   PetscErrorCode         ierr;
85*34a97f8cSStefano Zampini 
86*34a97f8cSStefano Zampini   PetscFunctionBegin;
87*34a97f8cSStefano Zampini   ierr = ISSorted(is_A_I,&is_sorted);CHKERRQ(ierr);
88*34a97f8cSStefano Zampini   if (!is_sorted) {
89*34a97f8cSStefano Zampini     SETERRQ(PetscObjectComm((PetscObject)is_A_I),PETSC_ERR_PLIB,"IS for I dofs should be shorted");
90*34a97f8cSStefano Zampini   }
91*34a97f8cSStefano Zampini   ierr = ISSorted(is_A_B,&is_sorted);CHKERRQ(ierr);
92*34a97f8cSStefano Zampini   if (!is_sorted) {
93*34a97f8cSStefano Zampini     SETERRQ(PetscObjectComm((PetscObject)is_A_B),PETSC_ERR_PLIB,"IS for B dofs should be shorted");
94*34a97f8cSStefano Zampini   }
95*34a97f8cSStefano Zampini 
96*34a97f8cSStefano Zampini   /* get sizes */
97*34a97f8cSStefano Zampini   ierr = ISGetLocalSize(is_A_I,&n_I);CHKERRQ(ierr);
98*34a97f8cSStefano Zampini   ierr = ISGetLocalSize(is_A_B,&n_B);CHKERRQ(ierr);
99*34a97f8cSStefano Zampini   n_local = n_I+n_B;
100*34a97f8cSStefano Zampini 
101*34a97f8cSStefano Zampini   /* maps */
102*34a97f8cSStefano Zampini   ierr = ISLocalToGlobalMappingCreateIS(is_A_B,&BtoNmap);CHKERRQ(ierr);
103*34a97f8cSStefano Zampini   if (nlayers >= 0 && xadj != NULL && adjncy != NULL) { /* I problems have a different size of the original ones */
104*34a97f8cSStefano Zampini     ierr = ISLocalToGlobalMappingCreateIS(is_A_I,&ItoNmap);CHKERRQ(ierr);
105*34a97f8cSStefano Zampini     /* allocate some auxiliary space */
106*34a97f8cSStefano Zampini     ierr = PetscMalloc1(n_local,&local_numbering);CHKERRQ(ierr);
107*34a97f8cSStefano Zampini     ierr = PetscBTCreate(n_local,&touched);CHKERRQ(ierr);
108*34a97f8cSStefano Zampini   } else {
109*34a97f8cSStefano Zampini     ItoNmap = 0;
110*34a97f8cSStefano Zampini     local_numbering = 0;
111*34a97f8cSStefano Zampini     touched = 0;
112*34a97f8cSStefano Zampini   }
113*34a97f8cSStefano Zampini 
114*34a97f8cSStefano Zampini   /* get Schur complement matrices */
115*34a97f8cSStefano Zampini   ierr = MatSchurComplementGetSubMatrices(S,&A_II,NULL,&A_IB,&A_BI,&A_BB);CHKERRQ(ierr);
116*34a97f8cSStefano Zampini 
117*34a97f8cSStefano Zampini   /* allocate space for schur complements */
118*34a97f8cSStefano Zampini   ierr = PetscMalloc5(ncc,&sub_schurs->is_AEj_I,ncc,&sub_schurs->is_AEj_B,ncc,&sub_schurs->S_Ej,ncc,&sub_schurs->work1,ncc,&sub_schurs->work2);CHKERRQ(ierr);
119*34a97f8cSStefano Zampini   sub_schurs->n_subs = ncc;
120*34a97f8cSStefano Zampini 
121*34a97f8cSStefano Zampini   /* cycle on subsets and extract schur complements */
122*34a97f8cSStefano Zampini   for (i=0;i<sub_schurs->n_subs;i++) {
123*34a97f8cSStefano Zampini     Mat      AE_II,AE_IE,AE_EI,AE_EE;
124*34a97f8cSStefano Zampini     IS       is_I,is_subset_B;
125*34a97f8cSStefano Zampini     PetscInt subset_size;
126*34a97f8cSStefano Zampini 
127*34a97f8cSStefano Zampini     /* size of subset */
128*34a97f8cSStefano Zampini     subset_size = cptr[i+1]-cptr[i];
129*34a97f8cSStefano Zampini     if (subset_size == 1 && skip_corners) {
130*34a97f8cSStefano Zampini       sub_schurs->is_AEj_B[i] = 0;
131*34a97f8cSStefano Zampini       sub_schurs->is_AEj_I[i] = 0;
132*34a97f8cSStefano Zampini       sub_schurs->S_Ej[i] = 0;
133*34a97f8cSStefano Zampini       sub_schurs->work1[i] = 0;
134*34a97f8cSStefano Zampini       sub_schurs->work2[i] = 0;
135*34a97f8cSStefano Zampini       continue;
136*34a97f8cSStefano Zampini     }
137*34a97f8cSStefano Zampini 
138*34a97f8cSStefano Zampini     /* get IS for subsets in original numbering and in B numbering */
139*34a97f8cSStefano Zampini     ierr = ISCreateGeneral(PetscObjectComm((PetscObject)BtoNmap),subset_size,queue+cptr[i],PETSC_COPY_VALUES,&sub_schurs->is_AEj_B[i]);CHKERRQ(ierr);
140*34a97f8cSStefano Zampini     ierr = ISSort(sub_schurs->is_AEj_B[i]);CHKERRQ(ierr);
141*34a97f8cSStefano Zampini     ierr = ISGlobalToLocalMappingApplyIS(BtoNmap,IS_GTOLM_DROP,sub_schurs->is_AEj_B[i],&is_subset_B);CHKERRQ(ierr);
142*34a97f8cSStefano Zampini 
143*34a97f8cSStefano Zampini     /* BB block on subset */
144*34a97f8cSStefano Zampini     ierr = MatGetSubMatrix(A_BB,is_subset_B,is_subset_B,MAT_INITIAL_MATRIX,&AE_EE);CHKERRQ(ierr);
145*34a97f8cSStefano Zampini 
146*34a97f8cSStefano Zampini     if (ItoNmap) { /* is ItoNmap has been computed, extracts only a part of I dofs */
147*34a97f8cSStefano Zampini       const PetscInt* idx_B;
148*34a97f8cSStefano Zampini       PetscInt        n_local_dofs,n_prev_added,j,layer;
149*34a97f8cSStefano Zampini 
150*34a97f8cSStefano Zampini       /* all boundary dofs must be skipped when adding layers */
151*34a97f8cSStefano Zampini       ierr = PetscBTMemzero(n_local,touched);CHKERRQ(ierr);
152*34a97f8cSStefano Zampini       ierr = ISGetIndices(is_A_B,&idx_B);CHKERRQ(ierr);
153*34a97f8cSStefano Zampini       for (j=0;j<n_B;j++) {
154*34a97f8cSStefano Zampini         ierr = PetscBTSet(touched,idx_B[j]);CHKERRQ(ierr);
155*34a97f8cSStefano Zampini       }
156*34a97f8cSStefano Zampini       ierr = ISRestoreIndices(is_A_B,&idx_B);CHKERRQ(ierr);
157*34a97f8cSStefano Zampini 
158*34a97f8cSStefano Zampini       /* add next layers of dofs */
159*34a97f8cSStefano Zampini       ierr = PetscMemcpy(local_numbering,queue+cptr[i],subset_size*sizeof(PetscInt));CHKERRQ(ierr);
160*34a97f8cSStefano Zampini       n_local_dofs = subset_size;
161*34a97f8cSStefano Zampini       n_prev_added = subset_size;
162*34a97f8cSStefano Zampini       for (layer=0;layer<nlayers;layer++) {
163*34a97f8cSStefano Zampini         PetscInt n_added;
164*34a97f8cSStefano Zampini         if (n_local_dofs == n_I+subset_size) break;
165*34a97f8cSStefano Zampini         if (n_local_dofs > n_I+subset_size) {
166*34a97f8cSStefano Zampini           SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error querying layer %d. Out of bound access (%d > %d)",layer,n_local_dofs,n_I+subset_size);
167*34a97f8cSStefano Zampini         }
168*34a97f8cSStefano Zampini         ierr = PCBDDCAdjGetNextLayer_Private(local_numbering+n_local_dofs,n_prev_added,touched,xadj,adjncy,&n_added);CHKERRQ(ierr);
169*34a97f8cSStefano Zampini         n_prev_added = n_added;
170*34a97f8cSStefano Zampini         n_local_dofs += n_added;
171*34a97f8cSStefano Zampini       }
172*34a97f8cSStefano Zampini 
173*34a97f8cSStefano Zampini       /* IS for I dofs in original numbering and in I numbering */
174*34a97f8cSStefano Zampini       ierr = ISCreateGeneral(PetscObjectComm((PetscObject)ItoNmap),n_local_dofs-subset_size,local_numbering+subset_size,PETSC_COPY_VALUES,&sub_schurs->is_AEj_I[i]);CHKERRQ(ierr);
175*34a97f8cSStefano Zampini       ierr = ISSort(sub_schurs->is_AEj_I[i]);CHKERRQ(ierr);
176*34a97f8cSStefano Zampini       ierr = ISGlobalToLocalMappingApplyIS(ItoNmap,IS_GTOLM_DROP,sub_schurs->is_AEj_I[i],&is_I);CHKERRQ(ierr);
177*34a97f8cSStefano Zampini 
178*34a97f8cSStefano Zampini       /* II block */
179*34a97f8cSStefano Zampini       ierr = MatGetSubMatrix(A_II,is_I,is_I,MAT_INITIAL_MATRIX,&AE_II);CHKERRQ(ierr);
180*34a97f8cSStefano Zampini     } else { /* in this case we can take references of already existing IS and matrices for I dofs */
181*34a97f8cSStefano Zampini       /* IS for I dofs in original numbering */
182*34a97f8cSStefano Zampini       ierr = PetscObjectReference((PetscObject)is_A_I);CHKERRQ(ierr);
183*34a97f8cSStefano Zampini       sub_schurs->is_AEj_I[i] = is_A_I;
184*34a97f8cSStefano Zampini 
185*34a97f8cSStefano Zampini       /* IS for I dofs in I numbering TODO: "first" argument of ISCreateStride is not general */
186*34a97f8cSStefano Zampini       ierr = ISCreateStride(PetscObjectComm((PetscObject)is_A_I),n_I,0,1,&is_I);CHKERRQ(ierr);
187*34a97f8cSStefano Zampini 
188*34a97f8cSStefano Zampini       /* II block is the same */
189*34a97f8cSStefano Zampini       ierr = PetscObjectReference((PetscObject)A_II);CHKERRQ(ierr);
190*34a97f8cSStefano Zampini       AE_II = A_II;
191*34a97f8cSStefano Zampini     }
192*34a97f8cSStefano Zampini 
193*34a97f8cSStefano Zampini     /* IE block */
194*34a97f8cSStefano Zampini     ierr = MatGetSubMatrix(A_IB,is_I,is_subset_B,MAT_INITIAL_MATRIX,&AE_IE);CHKERRQ(ierr);
195*34a97f8cSStefano Zampini 
196*34a97f8cSStefano Zampini     /* EI block */
197*34a97f8cSStefano Zampini     ierr = MatGetSubMatrix(A_BI,is_subset_B,is_I,MAT_INITIAL_MATRIX,&AE_EI);CHKERRQ(ierr);
198*34a97f8cSStefano Zampini 
199*34a97f8cSStefano Zampini     /* setup Schur complements on subset */
200*34a97f8cSStefano Zampini     ierr = MatCreateSchurComplement(AE_II,AE_II,AE_IE,AE_EI,AE_EE,&sub_schurs->S_Ej[i]);CHKERRQ(ierr);
201*34a97f8cSStefano Zampini     ierr = MatGetVecs(sub_schurs->S_Ej[i],&sub_schurs->work1[i],&sub_schurs->work2[i]);CHKERRQ(ierr);
202*34a97f8cSStefano Zampini     if (AE_II == A_II) { /* we can reuse the same ksp */
203*34a97f8cSStefano Zampini       KSP ksp;
204*34a97f8cSStefano Zampini       ierr = MatSchurComplementGetKSP(S,&ksp);CHKERRQ(ierr);
205*34a97f8cSStefano Zampini       ierr = MatSchurComplementSetKSP(sub_schurs->S_Ej[i],ksp);CHKERRQ(ierr);
206*34a97f8cSStefano Zampini     } else { /* build new ksp object which inherits ksp and pc types from the original one */
207*34a97f8cSStefano Zampini       KSP      origksp,schurksp;
208*34a97f8cSStefano Zampini       PC       origpc,schurpc;
209*34a97f8cSStefano Zampini       KSPType  ksp_type;
210*34a97f8cSStefano Zampini       PCType   pc_type;
211*34a97f8cSStefano Zampini       PetscInt n_internal;
212*34a97f8cSStefano Zampini 
213*34a97f8cSStefano Zampini       ierr = MatSchurComplementGetKSP(S,&origksp);CHKERRQ(ierr);
214*34a97f8cSStefano Zampini       ierr = MatSchurComplementGetKSP(sub_schurs->S_Ej[i],&schurksp);CHKERRQ(ierr);
215*34a97f8cSStefano Zampini       ierr = KSPGetType(origksp,&ksp_type);CHKERRQ(ierr);
216*34a97f8cSStefano Zampini       ierr = KSPSetType(schurksp,ksp_type);CHKERRQ(ierr);
217*34a97f8cSStefano Zampini       ierr = KSPGetPC(schurksp,&schurpc);CHKERRQ(ierr);
218*34a97f8cSStefano Zampini       ierr = KSPGetPC(origksp,&origpc);CHKERRQ(ierr);
219*34a97f8cSStefano Zampini       ierr = PCGetType(origpc,&pc_type);CHKERRQ(ierr);
220*34a97f8cSStefano Zampini       ierr = PCSetType(schurpc,pc_type);CHKERRQ(ierr);
221*34a97f8cSStefano Zampini       ierr = ISGetSize(is_I,&n_internal);CHKERRQ(ierr);
222*34a97f8cSStefano Zampini       if (n_internal) { /* UMFPACK gives error with 0 sized problems */
223*34a97f8cSStefano Zampini         MatSolverPackage solver=NULL;
224*34a97f8cSStefano Zampini         ierr = PCFactorGetMatSolverPackage(origpc,(const MatSolverPackage*)&solver);CHKERRQ(ierr);
225*34a97f8cSStefano Zampini         if (solver) {
226*34a97f8cSStefano Zampini           ierr = PCFactorSetMatSolverPackage(schurpc,solver);CHKERRQ(ierr);
227*34a97f8cSStefano Zampini         }
228*34a97f8cSStefano Zampini       }
229*34a97f8cSStefano Zampini       ierr = KSPSetUp(schurksp);CHKERRQ(ierr);
230*34a97f8cSStefano Zampini     }
231*34a97f8cSStefano Zampini     /* free */
232*34a97f8cSStefano Zampini     ierr = MatDestroy(&AE_II);CHKERRQ(ierr);
233*34a97f8cSStefano Zampini     ierr = MatDestroy(&AE_EE);CHKERRQ(ierr);
234*34a97f8cSStefano Zampini     ierr = MatDestroy(&AE_IE);CHKERRQ(ierr);
235*34a97f8cSStefano Zampini     ierr = MatDestroy(&AE_EI);CHKERRQ(ierr);
236*34a97f8cSStefano Zampini     ierr = ISDestroy(&is_I);CHKERRQ(ierr);
237*34a97f8cSStefano Zampini     ierr = ISDestroy(&is_subset_B);CHKERRQ(ierr);
238*34a97f8cSStefano Zampini   }
239*34a97f8cSStefano Zampini   /* free */
240*34a97f8cSStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&ItoNmap);CHKERRQ(ierr);
241*34a97f8cSStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&BtoNmap);CHKERRQ(ierr);
242*34a97f8cSStefano Zampini   ierr = PetscFree(local_numbering);CHKERRQ(ierr);
243*34a97f8cSStefano Zampini   ierr = PetscBTDestroy(&touched);CHKERRQ(ierr);
244*34a97f8cSStefano Zampini   PetscFunctionReturn(0);
245*34a97f8cSStefano Zampini }
246