1 #include <../src/ksp/pc/impls/bddc/bddc.h> 2 #include <../src/ksp/pc/impls/bddc/bddcprivate.h> 3 4 PETSC_STATIC_INLINE PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt*,PetscInt,PetscBT,PetscInt*,PetscInt*,PetscInt*); 5 static PetscErrorCode PCBDDCComputeExplicitSchur(Mat M, Mat *S); 6 7 #undef __FUNCT__ 8 #define __FUNCT__ "PCBDDCComputeExplicitSchur" 9 static PetscErrorCode PCBDDCComputeExplicitSchur(Mat M, Mat *S) 10 { 11 Mat B, C, D, Bd, Cd, AinvBd; 12 KSP ksp; 13 PC pc; 14 PetscBool isLU, isILU, isCHOL, Bdense, Cdense; 15 PetscReal fill = 2.0; 16 PetscMPIInt size; 17 PetscErrorCode ierr; 18 19 PetscFunctionBegin; 20 ierr = MPI_Comm_size(PetscObjectComm((PetscObject)M),&size);CHKERRQ(ierr); 21 if (size != 1) { 22 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not for parallel matrices"); 23 } 24 ierr = MatSchurComplementGetSubMatrices(M, NULL, NULL, &B, &C, &D);CHKERRQ(ierr); 25 ierr = MatSchurComplementGetKSP(M, &ksp);CHKERRQ(ierr); 26 ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr); 27 ierr = PetscObjectTypeCompare((PetscObject) pc, PCLU, &isLU);CHKERRQ(ierr); 28 ierr = PetscObjectTypeCompare((PetscObject) pc, PCILU, &isILU);CHKERRQ(ierr); 29 ierr = PetscObjectTypeCompare((PetscObject) pc, PCCHOLESKY, &isCHOL);CHKERRQ(ierr); 30 ierr = PetscObjectTypeCompare((PetscObject) B, MATSEQDENSE, &Bdense);CHKERRQ(ierr); 31 ierr = PetscObjectTypeCompare((PetscObject) C, MATSEQDENSE, &Cdense);CHKERRQ(ierr); 32 if (!Bdense) { 33 ierr = MatConvert(B, MATSEQDENSE, MAT_INITIAL_MATRIX, &Bd);CHKERRQ(ierr); 34 } else { 35 Bd = B; 36 } 37 38 if (isLU || isILU || isCHOL) { 39 Mat fact; 40 41 ierr = KSPSetUp(ksp);CHKERRQ(ierr); 42 ierr = PCFactorGetMatrix(pc, &fact);CHKERRQ(ierr); 43 ierr = MatDuplicate(Bd, MAT_DO_NOT_COPY_VALUES, &AinvBd);CHKERRQ(ierr); 44 ierr = MatMatSolve(fact, Bd, AinvBd);CHKERRQ(ierr); 45 } else { 46 Mat Ainvd; 47 48 ierr = PCComputeExplicitOperator(pc, &Ainvd);CHKERRQ(ierr); 49 ierr = MatMatMult(Ainvd, Bd, MAT_INITIAL_MATRIX, fill, &AinvBd);CHKERRQ(ierr); 50 ierr = MatDestroy(&Ainvd);CHKERRQ(ierr); 51 } 52 if (!Bdense) { 53 ierr = MatDestroy(&Bd);CHKERRQ(ierr); 54 } 55 if (!Cdense) { 56 ierr = MatConvert(C, MATSEQDENSE, MAT_INITIAL_MATRIX, &Cd);CHKERRQ(ierr); 57 } else { 58 Cd = C; 59 } 60 61 ierr = MatMatMult(Cd, AinvBd, MAT_INITIAL_MATRIX, fill, S);CHKERRQ(ierr); 62 ierr = MatDestroy(&AinvBd);CHKERRQ(ierr); 63 if (!Cdense) { 64 ierr = MatDestroy(&Cd);CHKERRQ(ierr); 65 } 66 67 if (D) { 68 Mat Dd; 69 PetscBool Ddense; 70 71 ierr = PetscObjectTypeCompare((PetscObject)D,MATSEQDENSE,&Ddense);CHKERRQ(ierr); 72 if (!Ddense) { 73 ierr = MatConvert(D, MATSEQDENSE, MAT_INITIAL_MATRIX, &Dd);CHKERRQ(ierr); 74 } else { 75 Dd = D; 76 } 77 ierr = MatAYPX(*S,-1.0,Dd,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 78 if (!Ddense) { 79 ierr = MatDestroy(&Dd);CHKERRQ(ierr); 80 } 81 } else { 82 ierr = MatScale(*S,-1.0);CHKERRQ(ierr); 83 } 84 PetscFunctionReturn(0); 85 } 86 87 #undef __FUNCT__ 88 #define __FUNCT__ "PCBDDCSubSchursSetUp" 89 PetscErrorCode PCBDDCSubSchursSetUp(PCBDDCSubSchurs sub_schurs, PetscInt xadj[], PetscInt adjncy[], PetscInt nlayers) 90 { 91 Mat A_II,A_IB,A_BI,A_BB; 92 Mat AE_II,*AE_IE,*AE_EI,*AE_EE; 93 Mat S_all,global_schur_subsets,work_mat; 94 ISLocalToGlobalMapping l2gmap_subsets; 95 IS is_I,*is_subset_B,temp_is; 96 PetscInt *nnz,*all_local_idx_G,*all_local_idx_B,*all_local_idx_N; 97 PetscInt i,subset_size,max_subset_size; 98 PetscInt extra,local_size,global_size; 99 PetscErrorCode ierr; 100 101 PetscFunctionBegin; 102 103 /* allocate space for schur complements */ 104 ierr = PetscMalloc4(sub_schurs->n_subs,&sub_schurs->is_AEj_B, 105 sub_schurs->n_subs,&sub_schurs->S_Ej, 106 sub_schurs->n_subs,&sub_schurs->work1, 107 sub_schurs->n_subs,&sub_schurs->work2);CHKERRQ(ierr); 108 109 /* get Schur complement matrices */ 110 if (!sub_schurs->use_mumps) { 111 ierr = MatSchurComplementGetSubMatrices(sub_schurs->S,&A_II,NULL,&A_IB,&A_BI,&A_BB);CHKERRQ(ierr); 112 ierr = PetscMalloc4(sub_schurs->n_subs,&is_subset_B, 113 sub_schurs->n_subs,&AE_IE, 114 sub_schurs->n_subs,&AE_EI, 115 sub_schurs->n_subs,&AE_EE);CHKERRQ(ierr); 116 } 117 118 /* determine interior problems */ 119 if (nlayers >= 0 && xadj != NULL && adjncy != NULL) { /* Interior problems can be different from the original one */ 120 PetscBT touched; 121 const PetscInt* idx_B; 122 PetscInt n_I,n_B,n_local_dofs,n_prev_added,j,layer,*local_numbering; 123 124 /* get sizes */ 125 ierr = ISGetLocalSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr); 126 ierr = ISGetLocalSize(sub_schurs->is_B,&n_B);CHKERRQ(ierr); 127 128 ierr = PetscMalloc1(n_I+n_B,&local_numbering);CHKERRQ(ierr); 129 ierr = PetscBTCreate(n_I+n_B,&touched);CHKERRQ(ierr); 130 ierr = PetscBTMemzero(n_I+n_B,touched);CHKERRQ(ierr); 131 132 /* all boundary dofs must be skipped when adding layers */ 133 ierr = ISGetIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr); 134 for (j=0;j<n_B;j++) { 135 ierr = PetscBTSet(touched,idx_B[j]);CHKERRQ(ierr); 136 } 137 ierr = PetscMemcpy(local_numbering,idx_B,n_B*sizeof(PetscInt));CHKERRQ(ierr); 138 ierr = ISRestoreIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr); 139 140 /* add prescribed number of layers of dofs */ 141 n_local_dofs = n_B; 142 n_prev_added = n_B; 143 for (layer=0;layer<nlayers;layer++) { 144 PetscInt n_added; 145 if (n_local_dofs == n_I+n_B) break; 146 if (n_local_dofs > n_I+n_B) { 147 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); 148 } 149 ierr = PCBDDCAdjGetNextLayer_Private(local_numbering+n_local_dofs,n_prev_added,touched,xadj,adjncy,&n_added);CHKERRQ(ierr); 150 n_prev_added = n_added; 151 n_local_dofs += n_added; 152 if (!n_added) break; 153 } 154 ierr = PetscBTDestroy(&touched);CHKERRQ(ierr); 155 156 /* IS for I layer dofs in original numbering */ 157 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); 158 ierr = PetscFree(local_numbering);CHKERRQ(ierr); 159 ierr = ISSort(sub_schurs->is_I_layer);CHKERRQ(ierr); 160 /* IS for I layer dofs in I numbering */ 161 if (!sub_schurs->use_mumps) { 162 ISLocalToGlobalMapping ItoNmap; 163 ierr = ISLocalToGlobalMappingCreateIS(sub_schurs->is_I,&ItoNmap);CHKERRQ(ierr); 164 ierr = ISGlobalToLocalMappingApplyIS(ItoNmap,IS_GTOLM_DROP,sub_schurs->is_I_layer,&is_I);CHKERRQ(ierr); 165 ierr = ISLocalToGlobalMappingDestroy(&ItoNmap);CHKERRQ(ierr); 166 167 /* II block */ 168 ierr = MatGetSubMatrix(A_II,is_I,is_I,MAT_INITIAL_MATRIX,&AE_II);CHKERRQ(ierr); 169 } 170 } else { 171 PetscInt n_I; 172 173 /* IS for I dofs in original numbering */ 174 ierr = PetscObjectReference((PetscObject)sub_schurs->is_I);CHKERRQ(ierr); 175 sub_schurs->is_I_layer = sub_schurs->is_I; 176 177 /* IS for I dofs in I numbering (strided 1) */ 178 if (!sub_schurs->use_mumps) { 179 ierr = ISGetSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr); 180 ierr = ISCreateStride(PetscObjectComm((PetscObject)sub_schurs->is_I),n_I,0,1,&is_I);CHKERRQ(ierr); 181 182 /* II block is the same */ 183 ierr = PetscObjectReference((PetscObject)A_II);CHKERRQ(ierr); 184 AE_II = A_II; 185 } 186 } 187 188 for (i=0;i<sub_schurs->n_subs;i++) { 189 ierr = ISDuplicate(sub_schurs->is_subs[i],&sub_schurs->is_AEj_B[i]);CHKERRQ(ierr); 190 ierr = ISSort(sub_schurs->is_AEj_B[i]);CHKERRQ(ierr); 191 } 192 193 /* Get info on subset sizes and sum of all subsets sizes */ 194 max_subset_size = 0; 195 local_size = 0; 196 for (i=0;i<sub_schurs->n_subs_seq;i++) { 197 PetscInt j = sub_schurs->index_sequential[i]; 198 ierr = ISGetLocalSize(sub_schurs->is_AEj_B[j],&subset_size);CHKERRQ(ierr); 199 max_subset_size = PetscMax(subset_size,max_subset_size); 200 local_size += subset_size; 201 } 202 203 /* Work arrays for local indices */ 204 ierr = PetscMalloc1(local_size,&all_local_idx_B);CHKERRQ(ierr); 205 extra = 0; 206 if (sub_schurs->use_mumps) { 207 ierr = ISGetLocalSize(sub_schurs->is_I_layer,&extra);CHKERRQ(ierr); 208 } 209 ierr = PetscMalloc1(local_size+extra,&all_local_idx_N);CHKERRQ(ierr); 210 if (extra) { 211 const PetscInt *idxs; 212 ierr = ISGetIndices(sub_schurs->is_I_layer,&idxs);CHKERRQ(ierr); 213 ierr = PetscMemcpy(all_local_idx_N,idxs,extra*sizeof(PetscInt));CHKERRQ(ierr); 214 ierr = ISRestoreIndices(sub_schurs->is_I_layer,&idxs);CHKERRQ(ierr); 215 } 216 ierr = PetscMalloc1(local_size,&nnz);CHKERRQ(ierr); 217 218 /* Get local indices in local numbering */ 219 local_size = 0; 220 for (i=0;i<sub_schurs->n_subs_seq;i++) { 221 PetscInt j; 222 const PetscInt *idxs; 223 224 PetscInt local_problem_index = sub_schurs->index_sequential[i]; 225 ierr = ISGetLocalSize(sub_schurs->is_AEj_B[local_problem_index],&subset_size);CHKERRQ(ierr); 226 ierr = ISGetIndices(sub_schurs->is_AEj_B[local_problem_index],&idxs);CHKERRQ(ierr); 227 /* subset indices in local numbering */ 228 ierr = PetscMemcpy(all_local_idx_N+local_size+extra,idxs,subset_size*sizeof(PetscInt));CHKERRQ(ierr); 229 ierr = ISRestoreIndices(sub_schurs->is_AEj_B[local_problem_index],&idxs);CHKERRQ(ierr); 230 for (j=0;j<subset_size;j++) nnz[local_size+j] = subset_size; 231 local_size += subset_size; 232 } 233 234 S_all = NULL; 235 if (!sub_schurs->use_mumps) { 236 /* subsets in original and boundary numbering */ 237 for (i=0;i<sub_schurs->n_subs;i++) { 238 ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,sub_schurs->is_AEj_B[i],&is_subset_B[i]);CHKERRQ(ierr); 239 } 240 241 /* EE block */ 242 for (i=0;i<sub_schurs->n_subs;i++) { 243 ierr = MatGetSubMatrix(A_BB,is_subset_B[i],is_subset_B[i],MAT_INITIAL_MATRIX,&AE_EE[i]);CHKERRQ(ierr); 244 } 245 /* IE block */ 246 for (i=0;i<sub_schurs->n_subs;i++) { 247 ierr = MatGetSubMatrix(A_IB,is_I,is_subset_B[i],MAT_INITIAL_MATRIX,&AE_IE[i]);CHKERRQ(ierr); 248 } 249 /* EI block */ 250 for (i=0;i<sub_schurs->n_subs;i++) { 251 ierr = MatGetSubMatrix(A_BI,is_subset_B[i],is_I,MAT_INITIAL_MATRIX,&AE_EI[i]);CHKERRQ(ierr); 252 } 253 254 /* setup Schur complements on subset */ 255 for (i=0;i<sub_schurs->n_subs;i++) { 256 ierr = MatCreateSchurComplement(AE_II,AE_II,AE_IE[i],AE_EI[i],AE_EE[i],&sub_schurs->S_Ej[i]);CHKERRQ(ierr); 257 if (AE_II == A_II) { /* we can reuse the same ksp */ 258 KSP ksp; 259 ierr = MatSchurComplementGetKSP(sub_schurs->S,&ksp);CHKERRQ(ierr); 260 ierr = MatSchurComplementSetKSP(sub_schurs->S_Ej[i],ksp);CHKERRQ(ierr); 261 } else { /* build new ksp object which inherits ksp and pc types from the original one */ 262 KSP origksp,schurksp; 263 PC origpc,schurpc; 264 KSPType ksp_type; 265 PCType pc_type; 266 PetscInt n_internal; 267 268 ierr = MatSchurComplementGetKSP(sub_schurs->S,&origksp);CHKERRQ(ierr); 269 ierr = MatSchurComplementGetKSP(sub_schurs->S_Ej[i],&schurksp);CHKERRQ(ierr); 270 ierr = KSPGetType(origksp,&ksp_type);CHKERRQ(ierr); 271 ierr = KSPSetType(schurksp,ksp_type);CHKERRQ(ierr); 272 ierr = KSPGetPC(schurksp,&schurpc);CHKERRQ(ierr); 273 ierr = KSPGetPC(origksp,&origpc);CHKERRQ(ierr); 274 ierr = PCGetType(origpc,&pc_type);CHKERRQ(ierr); 275 ierr = PCSetType(schurpc,pc_type);CHKERRQ(ierr); 276 ierr = ISGetSize(is_I,&n_internal);CHKERRQ(ierr); 277 if (n_internal) { /* UMFPACK gives error with 0 sized problems */ 278 MatSolverPackage solver=NULL; 279 ierr = PCFactorGetMatSolverPackage(origpc,(const MatSolverPackage*)&solver);CHKERRQ(ierr); 280 if (solver) { 281 ierr = PCFactorSetMatSolverPackage(schurpc,solver);CHKERRQ(ierr); 282 } 283 } 284 ierr = KSPSetUp(schurksp);CHKERRQ(ierr); 285 } 286 } 287 /* free */ 288 ierr = ISDestroy(&is_I);CHKERRQ(ierr); 289 ierr = MatDestroy(&AE_II);CHKERRQ(ierr); 290 for (i=0;i<sub_schurs->n_subs;i++) { 291 ierr = MatDestroy(&AE_EE[i]);CHKERRQ(ierr); 292 ierr = MatDestroy(&AE_IE[i]);CHKERRQ(ierr); 293 ierr = MatDestroy(&AE_EI[i]);CHKERRQ(ierr); 294 ierr = ISDestroy(&is_subset_B[i]);CHKERRQ(ierr); 295 } 296 ierr = PetscFree4(is_subset_B,AE_IE,AE_EI,AE_EE);CHKERRQ(ierr); 297 #if defined(PETSC_HAVE_MUMPS) 298 } else { 299 Mat A,F; 300 IS is_A_all; 301 PetscBool is_symmetric; 302 PetscInt *idxs_schur,n_I; 303 304 /* get working mat */ 305 ierr = ISGetLocalSize(sub_schurs->is_I_layer,&n_I);CHKERRQ(ierr); 306 ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size+n_I,all_local_idx_N,PETSC_COPY_VALUES,&is_A_all);CHKERRQ(ierr); 307 ierr = MatGetSubMatrixUnsorted(sub_schurs->A,is_A_all,is_A_all,&A);CHKERRQ(ierr); 308 ierr = ISDestroy(&is_A_all);CHKERRQ(ierr); 309 310 ierr = MatIsSymmetric(sub_schurs->A,0.0,&is_symmetric);CHKERRQ(ierr); 311 if (is_symmetric) { 312 ierr = MatGetFactor(A,MATSOLVERMUMPS,MAT_FACTOR_CHOLESKY,&F);CHKERRQ(ierr); 313 } else { 314 ierr = MatGetFactor(A,MATSOLVERMUMPS,MAT_FACTOR_LU,&F);CHKERRQ(ierr); 315 } 316 317 /* subsets ordered last */ 318 ierr = PetscMalloc1(local_size,&idxs_schur);CHKERRQ(ierr); 319 for (i=0;i<local_size;i++) { 320 idxs_schur[i] = n_I+i+1; 321 } 322 ierr = MatMumpsSetSchurIndices(F,local_size,idxs_schur);CHKERRQ(ierr); 323 ierr = PetscFree(idxs_schur);CHKERRQ(ierr); 324 325 /* factorization step */ 326 if (is_symmetric) { 327 ierr = MatCholeskyFactorSymbolic(F,A,NULL,NULL);CHKERRQ(ierr); 328 ierr = MatCholeskyFactorNumeric(F,A,NULL);CHKERRQ(ierr); 329 } else { 330 ierr = MatLUFactorSymbolic(F,A,NULL,NULL,NULL);CHKERRQ(ierr); 331 ierr = MatLUFactorNumeric(F,A,NULL);CHKERRQ(ierr); 332 } 333 334 /* get explicit Schur Complement computed during numeric factorization */ 335 ierr = MatMumpsGetSchurComplement(F,&S_all);CHKERRQ(ierr); 336 337 /* free workspace */ 338 ierr = MatDestroy(&A);CHKERRQ(ierr); 339 ierr = MatDestroy(&F);CHKERRQ(ierr); 340 341 /* unused */ 342 for (i=0;i<sub_schurs->n_subs;i++) { 343 sub_schurs->S_Ej[i] = 0; 344 } 345 #endif 346 } 347 348 /* TODO: just for compatibility with the previous version, needs to be fixed */ 349 for (i=0;i<sub_schurs->n_subs_par;i++) { 350 PetscInt j = sub_schurs->index_parallel[i]; 351 ierr = MatCreateVecs(sub_schurs->S_Ej[j],&sub_schurs->work1[j],&sub_schurs->work2[j]);CHKERRQ(ierr); 352 } 353 for (i=0;i<sub_schurs->n_subs_seq;i++) { 354 sub_schurs->work1[sub_schurs->index_sequential[i]] = 0; 355 sub_schurs->work2[sub_schurs->index_sequential[i]] = 0; 356 } 357 358 if (!sub_schurs->n_subs_seq_g) { 359 sub_schurs->S_Ej_all = 0; 360 sub_schurs->sum_S_Ej_all = 0; 361 sub_schurs->is_Ej_all = 0; 362 PetscFunctionReturn(0); 363 } 364 365 /* subset indices in local boundary numbering */ 366 ierr = ISGlobalToLocalMappingApply(sub_schurs->BtoNmap,IS_GTOLM_DROP,local_size,all_local_idx_N+extra,&subset_size,all_local_idx_B);CHKERRQ(ierr); 367 if (subset_size != local_size) { 368 SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in sub_schurs serial (BtoNmap)! %d != %d\n",subset_size,local_size); 369 } 370 ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size,all_local_idx_B,PETSC_OWN_POINTER,&sub_schurs->is_Ej_all);CHKERRQ(ierr); 371 372 /* Local matrix of all local Schur on subsets transposed */ 373 ierr = MatCreate(PETSC_COMM_SELF,&sub_schurs->S_Ej_all);CHKERRQ(ierr); 374 ierr = MatSetSizes(sub_schurs->S_Ej_all,PETSC_DECIDE,PETSC_DECIDE,local_size,local_size);CHKERRQ(ierr); 375 ierr = MatSetType(sub_schurs->S_Ej_all,MATAIJ);CHKERRQ(ierr); 376 ierr = MatSeqAIJSetPreallocation(sub_schurs->S_Ej_all,0,nnz);CHKERRQ(ierr); 377 ierr = PetscFree(nnz);CHKERRQ(ierr); 378 379 if (!sub_schurs->use_mumps) { 380 PetscScalar *fill_vals; 381 PetscInt *dummy_idx; 382 383 /* Work arrays */ 384 ierr = PetscMalloc1(max_subset_size,&dummy_idx);CHKERRQ(ierr); 385 386 /* Loop on local problems end compute Schur complements explicitly */ 387 local_size = 0; 388 for (i=0;i<sub_schurs->n_subs_seq;i++) { 389 Mat S_Ej_expl; 390 PetscInt j,lpi = sub_schurs->index_sequential[i]; 391 PetscBool Sdense; 392 393 ierr = PCBDDCComputeExplicitSchur(sub_schurs->S_Ej[lpi],&S_Ej_expl);CHKERRQ(ierr); 394 ierr = ISGetLocalSize(sub_schurs->is_AEj_B[lpi],&subset_size);CHKERRQ(ierr); 395 ierr = PetscObjectTypeCompare((PetscObject)S_Ej_expl,MATSEQDENSE,&Sdense);CHKERRQ(ierr); 396 if (Sdense) { 397 for (j=0;j<subset_size;j++) { 398 dummy_idx[j]=local_size+j; 399 } 400 ierr = MatDenseGetArray(S_Ej_expl,&fill_vals);CHKERRQ(ierr); 401 ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,fill_vals,INSERT_VALUES);CHKERRQ(ierr); 402 ierr = MatDenseRestoreArray(S_Ej_expl,&fill_vals);CHKERRQ(ierr); 403 } else { 404 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not yet implemented for sparse matrices"); 405 } 406 local_size += subset_size; 407 } 408 ierr = PetscFree(dummy_idx);CHKERRQ(ierr); 409 } else { 410 PetscInt *dummy_idx,n_all; 411 PetscScalar *vals,*fill_vals; 412 413 /* Work arrays */ 414 ierr = PetscMalloc1(max_subset_size,&dummy_idx);CHKERRQ(ierr); 415 ierr = MatGetSize(S_all,&n_all,NULL);CHKERRQ(ierr); 416 ierr = MatDenseGetArray(S_all,&vals);CHKERRQ(ierr); 417 local_size = 0; 418 subset_size = 0; 419 fill_vals = vals; 420 for (i=0;i<sub_schurs->n_subs_seq;i++) { 421 PetscInt j,lpi = sub_schurs->index_sequential[i]; 422 423 fill_vals += subset_size; 424 ierr = ISGetLocalSize(sub_schurs->is_AEj_B[lpi],&subset_size);CHKERRQ(ierr); 425 for (j=0;j<subset_size;j++) { 426 dummy_idx[j]=local_size+j; 427 } 428 for (j=0;j<subset_size;j++) { 429 ierr = MatSetValues(sub_schurs->S_Ej_all,1,dummy_idx+j,subset_size,dummy_idx,fill_vals,INSERT_VALUES);CHKERRQ(ierr); 430 fill_vals += n_all; 431 } 432 local_size += subset_size; 433 } 434 ierr = MatDenseRestoreArray(S_all,&vals);CHKERRQ(ierr); 435 ierr = PetscFree(dummy_idx);CHKERRQ(ierr); 436 } 437 ierr = MatAssemblyBegin(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 438 ierr = MatAssemblyEnd(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 439 440 /* Global matrix of all assembled Schur on subsets */ 441 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); 442 ierr = ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)sub_schurs->l2gmap),1,local_size,all_local_idx_G,PETSC_COPY_VALUES,&l2gmap_subsets);CHKERRQ(ierr); 443 ierr = MatCreateIS(PetscObjectComm((PetscObject)sub_schurs->l2gmap),1,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size,l2gmap_subsets,&work_mat);CHKERRQ(ierr); 444 ierr = ISLocalToGlobalMappingDestroy(&l2gmap_subsets);CHKERRQ(ierr); 445 ierr = MatISSetLocalMat(work_mat,sub_schurs->S_Ej_all);CHKERRQ(ierr); 446 ierr = MatISGetMPIXAIJ(work_mat,MAT_INITIAL_MATRIX,&global_schur_subsets);CHKERRQ(ierr); 447 448 /* Get local part of (\sum_j S_Ej) */ 449 ierr = PetscFree(all_local_idx_N);CHKERRQ(ierr); 450 ierr = ISCreateGeneral(PetscObjectComm((PetscObject)sub_schurs->l2gmap),local_size,all_local_idx_G,PETSC_OWN_POINTER,&temp_is);CHKERRQ(ierr); 451 ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr); 452 ierr = MatDestroy(&global_schur_subsets);CHKERRQ(ierr); 453 454 /* Computation of Stildas */ 455 if (S_all) { 456 PetscInt n_all; 457 PetscScalar *vals; 458 459 /* Work arrays */ 460 ierr = MatGetSize(S_all,&n_all,NULL);CHKERRQ(ierr); 461 ierr = MatDenseGetArray(S_all,&vals);CHKERRQ(ierr); 462 ierr = MatDenseRestoreArray(S_all,&vals);CHKERRQ(ierr); 463 } 464 465 ierr = MatDestroy(&S_all);CHKERRQ(ierr); 466 ierr = MatDestroy(&work_mat);CHKERRQ(ierr); 467 ierr = ISDestroy(&temp_is);CHKERRQ(ierr); 468 PetscFunctionReturn(0); 469 } 470 471 #undef __FUNCT__ 472 #define __FUNCT__ "PCBDDCSubSchursInit" 473 PetscErrorCode PCBDDCSubSchursInit(PCBDDCSubSchurs sub_schurs, Mat A, Mat S, IS is_I, IS is_B, PCBDDCGraph graph, ISLocalToGlobalMapping BtoNmap, PetscInt seqthreshold) 474 { 475 IS *faces,*edges,*all_cc; 476 PetscInt *index_sequential,*index_parallel; 477 PetscInt *auxlocal_sequential,*auxlocal_parallel; 478 PetscInt *auxglobal_sequential,*auxglobal_parallel; 479 PetscInt *auxmapping;//,*idxs; 480 PetscInt i,max_subset_size; 481 PetscInt n_sequential_problems,n_local_sequential_problems,n_parallel_problems,n_local_parallel_problems; 482 PetscInt n_faces,n_edges,n_all_cc; 483 PetscBool is_sorted; 484 PetscErrorCode ierr; 485 486 PetscFunctionBegin; 487 ierr = ISSorted(is_I,&is_sorted);CHKERRQ(ierr); 488 if (!is_sorted) { 489 SETERRQ(PetscObjectComm((PetscObject)is_I),PETSC_ERR_PLIB,"IS for I dofs should be shorted"); 490 } 491 ierr = ISSorted(is_B,&is_sorted);CHKERRQ(ierr); 492 if (!is_sorted) { 493 SETERRQ(PetscObjectComm((PetscObject)is_B),PETSC_ERR_PLIB,"IS for B dofs should be shorted"); 494 } 495 496 /* reset any previous data */ 497 ierr = PCBDDCSubSchursReset(sub_schurs);CHKERRQ(ierr); 498 499 /* get index sets for faces and edges */ 500 ierr = PCBDDCGraphGetCandidatesIS(graph,&n_faces,&faces,&n_edges,&edges,NULL);CHKERRQ(ierr); 501 n_all_cc = n_faces+n_edges; 502 ierr = PetscMalloc1(n_all_cc,&all_cc);CHKERRQ(ierr); 503 for (i=0;i<n_faces;i++) { 504 all_cc[i] = faces[i]; 505 } 506 for (i=0;i<n_edges;i++) { 507 all_cc[n_faces+i] = edges[i]; 508 } 509 ierr = PetscFree(faces);CHKERRQ(ierr); 510 ierr = PetscFree(edges);CHKERRQ(ierr); 511 512 /* map interface's subsets */ 513 max_subset_size = 0; 514 for (i=0;i<n_all_cc;i++) { 515 PetscInt subset_size; 516 ierr = ISGetLocalSize(all_cc[i],&subset_size);CHKERRQ(ierr); 517 max_subset_size = PetscMax(max_subset_size,subset_size); 518 } 519 ierr = PetscMalloc1(max_subset_size,&auxmapping);CHKERRQ(ierr); 520 ierr = PetscMalloc2(graph->ncc,&auxlocal_sequential, 521 graph->ncc,&auxlocal_parallel);CHKERRQ(ierr); 522 ierr = PetscMalloc2(graph->ncc,&index_sequential, 523 graph->ncc,&index_parallel);CHKERRQ(ierr); 524 525 /* if threshold is negative uses all sequential problems (possibly using MUMPS) */ 526 sub_schurs->use_mumps = PETSC_FALSE; 527 if (seqthreshold < 0) { 528 seqthreshold = max_subset_size; 529 #if defined(PETSC_HAVE_MUMPS) 530 sub_schurs->use_mumps = !!A; 531 #endif 532 } 533 534 535 /* determine which problem has to be solved in parallel or sequentially */ 536 n_local_sequential_problems = 0; 537 n_local_parallel_problems = 0; 538 for (i=0;i<n_all_cc;i++) { 539 PetscInt subset_size,j,min_loc = 0; 540 const PetscInt *idxs; 541 542 ierr = ISGetLocalSize(all_cc[i],&subset_size);CHKERRQ(ierr); 543 ierr = ISGetIndices(all_cc[i],&idxs);CHKERRQ(ierr); 544 ierr = ISLocalToGlobalMappingApply(graph->l2gmap,subset_size,idxs,auxmapping);CHKERRQ(ierr); 545 for (j=1;j<subset_size;j++) { 546 if (auxmapping[j]<auxmapping[min_loc]) { 547 min_loc = j; 548 } 549 } 550 if (subset_size > seqthreshold) { 551 index_parallel[n_local_parallel_problems] = i; 552 auxlocal_parallel[n_local_parallel_problems] = idxs[min_loc]; 553 n_local_parallel_problems++; 554 } else { 555 index_sequential[n_local_sequential_problems] = i; 556 auxlocal_sequential[n_local_sequential_problems] = idxs[min_loc]; 557 n_local_sequential_problems++; 558 } 559 ierr = ISRestoreIndices(all_cc[i],&idxs);CHKERRQ(ierr); 560 } 561 562 /* Number parallel problems */ 563 auxglobal_parallel = 0; 564 ierr = PCBDDCSubsetNumbering(PetscObjectComm((PetscObject)graph->l2gmap),graph->l2gmap,n_local_parallel_problems,auxlocal_parallel,PETSC_NULL,&n_parallel_problems,&auxglobal_parallel);CHKERRQ(ierr); 565 566 /* Number sequential problems */ 567 auxglobal_sequential = 0; 568 ierr = PCBDDCSubsetNumbering(PetscObjectComm((PetscObject)graph->l2gmap),graph->l2gmap,n_local_sequential_problems,auxlocal_sequential,PETSC_NULL,&n_sequential_problems,&auxglobal_sequential);CHKERRQ(ierr); 569 570 /* update info in sub_schurs */ 571 if (sub_schurs->use_mumps && A) { 572 ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 573 sub_schurs->A = A; 574 } 575 ierr = PetscObjectReference((PetscObject)S);CHKERRQ(ierr); 576 sub_schurs->S = S; 577 ierr = PetscObjectReference((PetscObject)is_I);CHKERRQ(ierr); 578 sub_schurs->is_I = is_I; 579 ierr = PetscObjectReference((PetscObject)is_B);CHKERRQ(ierr); 580 sub_schurs->is_B = is_B; 581 ierr = PetscObjectReference((PetscObject)graph->l2gmap);CHKERRQ(ierr); 582 sub_schurs->l2gmap = graph->l2gmap; 583 ierr = PetscObjectReference((PetscObject)BtoNmap);CHKERRQ(ierr); 584 sub_schurs->BtoNmap = BtoNmap; 585 sub_schurs->n_subs_seq = n_local_sequential_problems; 586 sub_schurs->n_subs_par = n_local_parallel_problems; 587 sub_schurs->n_subs_seq_g = n_sequential_problems; 588 sub_schurs->n_subs_par_g = n_parallel_problems; 589 sub_schurs->n_subs = sub_schurs->n_subs_seq + sub_schurs->n_subs_par; 590 sub_schurs->is_subs = all_cc; 591 sub_schurs->index_sequential = index_sequential; 592 sub_schurs->index_parallel = index_parallel; 593 sub_schurs->auxglobal_sequential = auxglobal_sequential; 594 sub_schurs->auxglobal_parallel = auxglobal_parallel; 595 596 /* free workspace */ 597 ierr = PetscFree(auxmapping);CHKERRQ(ierr); 598 ierr = PetscFree2(auxlocal_sequential,auxlocal_parallel);CHKERRQ(ierr); 599 600 PetscFunctionReturn(0); 601 } 602 603 #undef __FUNCT__ 604 #define __FUNCT__ "PCBDDCSubSchursCreate" 605 PetscErrorCode PCBDDCSubSchursCreate(PCBDDCSubSchurs *sub_schurs) 606 { 607 PCBDDCSubSchurs schurs_ctx; 608 PetscErrorCode ierr; 609 610 PetscFunctionBegin; 611 ierr = PetscNew(&schurs_ctx);CHKERRQ(ierr); 612 schurs_ctx->n_subs = 0; 613 *sub_schurs = schurs_ctx; 614 PetscFunctionReturn(0); 615 } 616 617 #undef __FUNCT__ 618 #define __FUNCT__ "PCBDDCSubSchursDestroy" 619 PetscErrorCode PCBDDCSubSchursDestroy(PCBDDCSubSchurs *sub_schurs) 620 { 621 PetscErrorCode ierr; 622 623 PetscFunctionBegin; 624 ierr = PCBDDCSubSchursReset(*sub_schurs);CHKERRQ(ierr); 625 ierr = PetscFree(*sub_schurs);CHKERRQ(ierr); 626 PetscFunctionReturn(0); 627 } 628 629 #undef __FUNCT__ 630 #define __FUNCT__ "PCBDDCSubSchursReset" 631 PetscErrorCode PCBDDCSubSchursReset(PCBDDCSubSchurs sub_schurs) 632 { 633 PetscInt i; 634 PetscErrorCode ierr; 635 636 PetscFunctionBegin; 637 ierr = MatDestroy(&sub_schurs->A);CHKERRQ(ierr); 638 ierr = MatDestroy(&sub_schurs->S);CHKERRQ(ierr); 639 ierr = ISDestroy(&sub_schurs->is_I);CHKERRQ(ierr); 640 ierr = ISDestroy(&sub_schurs->is_B);CHKERRQ(ierr); 641 ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->l2gmap);CHKERRQ(ierr); 642 ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->BtoNmap);CHKERRQ(ierr); 643 ierr = MatDestroy(&sub_schurs->S_Ej_all);CHKERRQ(ierr); 644 ierr = MatDestroy(&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr); 645 ierr = ISDestroy(&sub_schurs->is_Ej_all);CHKERRQ(ierr); 646 for (i=0;i<sub_schurs->n_subs;i++) { 647 ierr = ISDestroy(&sub_schurs->is_subs[i]);CHKERRQ(ierr); 648 ierr = ISDestroy(&sub_schurs->is_AEj_B[i]);CHKERRQ(ierr); 649 ierr = MatDestroy(&sub_schurs->S_Ej[i]);CHKERRQ(ierr); 650 ierr = VecDestroy(&sub_schurs->work1[i]);CHKERRQ(ierr); 651 ierr = VecDestroy(&sub_schurs->work2[i]);CHKERRQ(ierr); 652 } 653 ierr = ISDestroy(&sub_schurs->is_I_layer);CHKERRQ(ierr); 654 if (sub_schurs->n_subs) { 655 ierr = PetscFree(sub_schurs->is_subs);CHKERRQ(ierr); 656 ierr = PetscFree4(sub_schurs->is_AEj_B,sub_schurs->S_Ej,sub_schurs->work1,sub_schurs->work2);CHKERRQ(ierr); 657 ierr = PetscFree2(sub_schurs->index_sequential,sub_schurs->index_parallel);CHKERRQ(ierr); 658 ierr = PetscFree(sub_schurs->auxglobal_sequential);CHKERRQ(ierr); 659 ierr = PetscFree(sub_schurs->auxglobal_parallel);CHKERRQ(ierr); 660 } 661 sub_schurs->n_subs = 0; 662 PetscFunctionReturn(0); 663 } 664 665 #undef __FUNCT__ 666 #define __FUNCT__ "PCBDDCAdjGetNextLayer_Private" 667 PETSC_STATIC_INLINE PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt* queue_tip,PetscInt n_prev,PetscBT touched,PetscInt* xadj,PetscInt* adjncy,PetscInt* n_added) 668 { 669 PetscInt i,j,n; 670 PetscErrorCode ierr; 671 672 PetscFunctionBegin; 673 n = 0; 674 for (i=-n_prev;i<0;i++) { 675 PetscInt start_dof = queue_tip[i]; 676 for (j=xadj[start_dof];j<xadj[start_dof+1];j++) { 677 PetscInt dof = adjncy[j]; 678 if (!PetscBTLookup(touched,dof)) { 679 ierr = PetscBTSet(touched,dof);CHKERRQ(ierr); 680 queue_tip[n] = dof; 681 n++; 682 } 683 } 684 } 685 *n_added = n; 686 PetscFunctionReturn(0); 687 } 688