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