1 #include <../src/ksp/pc/impls/bddc/bddc.h> 2 #include <../src/ksp/pc/impls/bddc/bddcprivate.h> 3 #include <petscblaslapack.h> 4 5 PETSC_STATIC_INLINE PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt*,PetscInt,PetscBT,PetscInt*,PetscInt*,PetscInt*); 6 static PetscErrorCode PCBDDCComputeExplicitSchur(Mat M, Mat *S); 7 8 #undef __FUNCT__ 9 #define __FUNCT__ "PCBDDCComputeExplicitSchur" 10 static PetscErrorCode PCBDDCComputeExplicitSchur(Mat M, Mat *S) 11 { 12 Mat B, C, D, Bd, Cd, AinvBd; 13 KSP ksp; 14 PC pc; 15 PetscBool isLU, isILU, isCHOL, Bdense, Cdense; 16 PetscReal fill = 2.0; 17 PetscMPIInt size; 18 PetscErrorCode ierr; 19 20 PetscFunctionBegin; 21 ierr = MPI_Comm_size(PetscObjectComm((PetscObject)M),&size);CHKERRQ(ierr); 22 if (size != 1) { 23 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not for parallel matrices"); 24 } 25 ierr = MatSchurComplementGetSubMatrices(M, NULL, NULL, &B, &C, &D);CHKERRQ(ierr); 26 ierr = MatSchurComplementGetKSP(M, &ksp);CHKERRQ(ierr); 27 ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr); 28 ierr = PetscObjectTypeCompare((PetscObject) pc, PCLU, &isLU);CHKERRQ(ierr); 29 ierr = PetscObjectTypeCompare((PetscObject) pc, PCILU, &isILU);CHKERRQ(ierr); 30 ierr = PetscObjectTypeCompare((PetscObject) pc, PCCHOLESKY, &isCHOL);CHKERRQ(ierr); 31 ierr = PetscObjectTypeCompare((PetscObject) B, MATSEQDENSE, &Bdense);CHKERRQ(ierr); 32 ierr = PetscObjectTypeCompare((PetscObject) C, MATSEQDENSE, &Cdense);CHKERRQ(ierr); 33 if (!Bdense) { 34 ierr = MatConvert(B, MATSEQDENSE, MAT_INITIAL_MATRIX, &Bd);CHKERRQ(ierr); 35 } else { 36 Bd = B; 37 } 38 39 if (isLU || isILU || isCHOL) { 40 Mat fact; 41 42 ierr = KSPSetUp(ksp);CHKERRQ(ierr); 43 ierr = PCFactorGetMatrix(pc, &fact);CHKERRQ(ierr); 44 ierr = MatDuplicate(Bd, MAT_DO_NOT_COPY_VALUES, &AinvBd);CHKERRQ(ierr); 45 ierr = MatMatSolve(fact, Bd, AinvBd);CHKERRQ(ierr); 46 } else { 47 Mat Ainvd; 48 49 ierr = PCComputeExplicitOperator(pc, &Ainvd);CHKERRQ(ierr); 50 ierr = MatMatMult(Ainvd, Bd, MAT_INITIAL_MATRIX, fill, &AinvBd);CHKERRQ(ierr); 51 ierr = MatDestroy(&Ainvd);CHKERRQ(ierr); 52 } 53 if (!Bdense) { 54 ierr = MatDestroy(&Bd);CHKERRQ(ierr); 55 } 56 if (!Cdense) { 57 ierr = MatConvert(C, MATSEQDENSE, MAT_INITIAL_MATRIX, &Cd);CHKERRQ(ierr); 58 } else { 59 Cd = C; 60 } 61 62 ierr = MatMatMult(Cd, AinvBd, MAT_INITIAL_MATRIX, fill, S);CHKERRQ(ierr); 63 ierr = MatDestroy(&AinvBd);CHKERRQ(ierr); 64 if (!Cdense) { 65 ierr = MatDestroy(&Cd);CHKERRQ(ierr); 66 } 67 68 if (D) { 69 Mat Dd; 70 PetscBool Ddense; 71 72 ierr = PetscObjectTypeCompare((PetscObject)D,MATSEQDENSE,&Ddense);CHKERRQ(ierr); 73 if (!Ddense) { 74 ierr = MatConvert(D, MATSEQDENSE, MAT_INITIAL_MATRIX, &Dd);CHKERRQ(ierr); 75 } else { 76 Dd = D; 77 } 78 ierr = MatAYPX(*S,-1.0,Dd,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 79 if (!Ddense) { 80 ierr = MatDestroy(&Dd);CHKERRQ(ierr); 81 } 82 } else { 83 ierr = MatScale(*S,-1.0);CHKERRQ(ierr); 84 } 85 PetscFunctionReturn(0); 86 } 87 88 #undef __FUNCT__ 89 #define __FUNCT__ "PCBDDCSubSchursSetUp" 90 PetscErrorCode PCBDDCSubSchursSetUp(PCBDDCSubSchurs sub_schurs, PetscInt xadj[], PetscInt adjncy[], PetscInt nlayers, PetscBool compute_Stilda, PetscBool deluxe, PetscBool use_edges, PetscBool use_faces) 91 { 92 Mat A_II,A_IB,A_BI,A_BB; 93 Mat AE_II,*AE_IE,*AE_EI,*AE_EE; 94 Mat S_all,global_schur_subsets,work_mat; 95 Mat S_Ej_tilda_all,S_Ej_inv_all; 96 ISLocalToGlobalMapping l2gmap_subsets; 97 IS is_I,*is_subset_B,temp_is; 98 PetscInt *nnz,*all_local_idx_G,*all_local_idx_B,*all_local_idx_N; 99 PetscInt i,subset_size,max_subset_size; 100 PetscInt extra,local_size,global_size; 101 PetscBool is_symmetric,Stilda_computed; 102 PetscBLASInt B_N,B_ierr,B_lwork,*pivots; 103 PetscScalar *work; 104 PetscErrorCode ierr; 105 106 PetscFunctionBegin; 107 /* get Schur complement matrices */ 108 if (!sub_schurs->use_mumps) { 109 if (compute_Stilda) { 110 SETERRQ(PetscObjectComm((PetscObject)sub_schurs->l2gmap),PETSC_ERR_SUP,"Computation of Stildas requires MUMPS"); 111 } 112 ierr = MatSchurComplementGetSubMatrices(sub_schurs->S,&A_II,NULL,&A_IB,&A_BI,&A_BB);CHKERRQ(ierr); 113 ierr = PetscMalloc4(sub_schurs->n_subs,&is_subset_B, 114 sub_schurs->n_subs,&AE_IE, 115 sub_schurs->n_subs,&AE_EI, 116 sub_schurs->n_subs,&AE_EE);CHKERRQ(ierr); 117 } else { 118 is_subset_B = NULL; 119 AE_IE = NULL; 120 AE_EI = NULL; 121 AE_EE = NULL; 122 } 123 124 /* determine interior problems */ 125 ierr = ISDestroy(&sub_schurs->is_I_layer);CHKERRQ(ierr); 126 if (nlayers >= 0 && xadj != NULL && adjncy != NULL) { /* Interior problems can be different from the original one */ 127 PetscBT touched; 128 const PetscInt* idx_B; 129 PetscInt n_I,n_B,n_local_dofs,n_prev_added,j,layer,*local_numbering; 130 131 /* get sizes */ 132 ierr = ISGetLocalSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr); 133 ierr = ISGetLocalSize(sub_schurs->is_B,&n_B);CHKERRQ(ierr); 134 135 ierr = PetscMalloc1(n_I+n_B,&local_numbering);CHKERRQ(ierr); 136 ierr = PetscBTCreate(n_I+n_B,&touched);CHKERRQ(ierr); 137 ierr = PetscBTMemzero(n_I+n_B,touched);CHKERRQ(ierr); 138 139 /* all boundary dofs must be skipped when adding layers */ 140 ierr = ISGetIndices(sub_schurs->is_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 = PetscMemcpy(local_numbering,idx_B,n_B*sizeof(PetscInt));CHKERRQ(ierr); 145 ierr = ISRestoreIndices(sub_schurs->is_B,&idx_B);CHKERRQ(ierr); 146 147 /* add prescribed number of layers of dofs */ 148 n_local_dofs = n_B; 149 n_prev_added = n_B; 150 for (layer=0;layer<nlayers;layer++) { 151 PetscInt n_added; 152 if (n_local_dofs == n_I+n_B) break; 153 if (n_local_dofs > n_I+n_B) { 154 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); 155 } 156 ierr = PCBDDCAdjGetNextLayer_Private(local_numbering+n_local_dofs,n_prev_added,touched,xadj,adjncy,&n_added);CHKERRQ(ierr); 157 n_prev_added = n_added; 158 n_local_dofs += n_added; 159 if (!n_added) break; 160 } 161 ierr = PetscBTDestroy(&touched);CHKERRQ(ierr); 162 163 /* IS for I layer dofs in original numbering */ 164 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); 165 ierr = PetscFree(local_numbering);CHKERRQ(ierr); 166 ierr = ISSort(sub_schurs->is_I_layer);CHKERRQ(ierr); 167 /* IS for I layer dofs in I numbering */ 168 if (!sub_schurs->use_mumps) { 169 ISLocalToGlobalMapping ItoNmap; 170 ierr = ISLocalToGlobalMappingCreateIS(sub_schurs->is_I,&ItoNmap);CHKERRQ(ierr); 171 ierr = ISGlobalToLocalMappingApplyIS(ItoNmap,IS_GTOLM_DROP,sub_schurs->is_I_layer,&is_I);CHKERRQ(ierr); 172 ierr = ISLocalToGlobalMappingDestroy(&ItoNmap);CHKERRQ(ierr); 173 174 /* II block */ 175 ierr = MatGetSubMatrix(A_II,is_I,is_I,MAT_INITIAL_MATRIX,&AE_II);CHKERRQ(ierr); 176 } 177 } else { 178 PetscInt n_I; 179 180 /* IS for I dofs in original numbering */ 181 ierr = PetscObjectReference((PetscObject)sub_schurs->is_I);CHKERRQ(ierr); 182 sub_schurs->is_I_layer = sub_schurs->is_I; 183 184 /* IS for I dofs in I numbering (strided 1) */ 185 if (!sub_schurs->use_mumps) { 186 ierr = ISGetSize(sub_schurs->is_I,&n_I);CHKERRQ(ierr); 187 ierr = ISCreateStride(PetscObjectComm((PetscObject)sub_schurs->is_I),n_I,0,1,&is_I);CHKERRQ(ierr); 188 189 /* II block is the same */ 190 ierr = PetscObjectReference((PetscObject)A_II);CHKERRQ(ierr); 191 AE_II = A_II; 192 } 193 } 194 195 /* Get info on subset sizes and sum of all subsets sizes */ 196 max_subset_size = 0; 197 local_size = 0; 198 for (i=0;i<sub_schurs->n_subs_seq;i++) { 199 PetscInt j = sub_schurs->index_sequential[i]; 200 ierr = ISGetLocalSize(sub_schurs->is_subs[j],&subset_size);CHKERRQ(ierr); 201 max_subset_size = PetscMax(subset_size,max_subset_size); 202 local_size += subset_size; 203 } 204 205 /* Work arrays for local indices */ 206 ierr = PetscMalloc1(local_size,&all_local_idx_B);CHKERRQ(ierr); 207 extra = 0; 208 if (sub_schurs->use_mumps) { 209 ierr = ISGetLocalSize(sub_schurs->is_I_layer,&extra);CHKERRQ(ierr); 210 } 211 ierr = PetscMalloc1(local_size+extra,&all_local_idx_N);CHKERRQ(ierr); 212 if (extra) { 213 const PetscInt *idxs; 214 ierr = ISGetIndices(sub_schurs->is_I_layer,&idxs);CHKERRQ(ierr); 215 ierr = PetscMemcpy(all_local_idx_N,idxs,extra*sizeof(PetscInt));CHKERRQ(ierr); 216 ierr = ISRestoreIndices(sub_schurs->is_I_layer,&idxs);CHKERRQ(ierr); 217 } 218 ierr = PetscMalloc1(local_size,&nnz);CHKERRQ(ierr); 219 220 /* Get local indices in local numbering */ 221 local_size = 0; 222 for (i=0;i<sub_schurs->n_subs_seq;i++) { 223 PetscInt j; 224 const PetscInt *idxs; 225 226 PetscInt local_problem_index = sub_schurs->index_sequential[i]; 227 ierr = ISGetLocalSize(sub_schurs->is_subs[local_problem_index],&subset_size);CHKERRQ(ierr); 228 ierr = ISGetIndices(sub_schurs->is_subs[local_problem_index],&idxs);CHKERRQ(ierr); 229 /* subset indices in local numbering */ 230 ierr = PetscMemcpy(all_local_idx_N+local_size+extra,idxs,subset_size*sizeof(PetscInt));CHKERRQ(ierr); 231 ierr = ISRestoreIndices(sub_schurs->is_subs[local_problem_index],&idxs);CHKERRQ(ierr); 232 for (j=0;j<subset_size;j++) nnz[local_size+j] = subset_size; 233 local_size += subset_size; 234 } 235 236 S_all = NULL; 237 ierr = MatIsSymmetric(sub_schurs->A,0.0,&is_symmetric);CHKERRQ(ierr); 238 if (!sub_schurs->use_mumps) { 239 /* subsets in original and boundary numbering */ 240 for (i=0;i<sub_schurs->n_subs;i++) { 241 ierr = ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap,IS_GTOLM_DROP,sub_schurs->is_subs[i],&is_subset_B[i]);CHKERRQ(ierr); 242 } 243 244 /* EE block */ 245 for (i=0;i<sub_schurs->n_subs;i++) { 246 ierr = MatGetSubMatrix(A_BB,is_subset_B[i],is_subset_B[i],MAT_INITIAL_MATRIX,&AE_EE[i]);CHKERRQ(ierr); 247 } 248 /* IE block */ 249 for (i=0;i<sub_schurs->n_subs;i++) { 250 ierr = MatGetSubMatrix(A_IB,is_I,is_subset_B[i],MAT_INITIAL_MATRIX,&AE_IE[i]);CHKERRQ(ierr); 251 } 252 /* EI block */ 253 for (i=0;i<sub_schurs->n_subs;i++) { 254 ierr = MatGetSubMatrix(A_BI,is_subset_B[i],is_I,MAT_INITIAL_MATRIX,&AE_EI[i]);CHKERRQ(ierr); 255 } 256 257 /* setup Schur complements on subset */ 258 for (i=0;i<sub_schurs->n_subs;i++) { 259 ierr = MatDestroy(&sub_schurs->S_Ej[i]);CHKERRQ(ierr); 260 ierr = MatCreateSchurComplement(AE_II,AE_II,AE_IE[i],AE_EI[i],AE_EE[i],&sub_schurs->S_Ej[i]);CHKERRQ(ierr); 261 ierr = MatDestroy(&AE_EE[i]);CHKERRQ(ierr); 262 ierr = MatDestroy(&AE_IE[i]);CHKERRQ(ierr); 263 ierr = MatDestroy(&AE_EI[i]);CHKERRQ(ierr); 264 if (AE_II == A_II) { /* we can reuse the same ksp */ 265 KSP ksp; 266 ierr = MatSchurComplementGetKSP(sub_schurs->S,&ksp);CHKERRQ(ierr); 267 ierr = MatSchurComplementSetKSP(sub_schurs->S_Ej[i],ksp);CHKERRQ(ierr); 268 } else { /* build new ksp object which inherits ksp and pc types from the original one */ 269 KSP origksp,schurksp; 270 PC origpc,schurpc; 271 KSPType ksp_type; 272 PCType pc_type; 273 PetscInt n_internal; 274 275 ierr = MatSchurComplementGetKSP(sub_schurs->S,&origksp);CHKERRQ(ierr); 276 ierr = MatSchurComplementGetKSP(sub_schurs->S_Ej[i],&schurksp);CHKERRQ(ierr); 277 ierr = KSPGetType(origksp,&ksp_type);CHKERRQ(ierr); 278 ierr = KSPSetType(schurksp,ksp_type);CHKERRQ(ierr); 279 ierr = KSPGetPC(schurksp,&schurpc);CHKERRQ(ierr); 280 ierr = KSPGetPC(origksp,&origpc);CHKERRQ(ierr); 281 ierr = PCGetType(origpc,&pc_type);CHKERRQ(ierr); 282 ierr = PCSetType(schurpc,pc_type);CHKERRQ(ierr); 283 ierr = ISGetSize(is_I,&n_internal);CHKERRQ(ierr); 284 if (n_internal) { /* UMFPACK gives error with 0 sized problems */ 285 MatSolverPackage solver=NULL; 286 ierr = PCFactorGetMatSolverPackage(origpc,(const MatSolverPackage*)&solver);CHKERRQ(ierr); 287 if (solver) { 288 ierr = PCFactorSetMatSolverPackage(schurpc,solver);CHKERRQ(ierr); 289 } 290 } 291 ierr = KSPSetUp(schurksp);CHKERRQ(ierr); 292 } 293 } 294 /* free */ 295 ierr = ISDestroy(&is_I);CHKERRQ(ierr); 296 ierr = MatDestroy(&AE_II);CHKERRQ(ierr); 297 for (i=0;i<sub_schurs->n_subs;i++) { 298 ierr = ISDestroy(&is_subset_B[i]);CHKERRQ(ierr); 299 } 300 ierr = PetscFree4(is_subset_B,AE_IE,AE_EI,AE_EE);CHKERRQ(ierr); 301 #if defined(PETSC_HAVE_MUMPS) 302 } else { 303 Mat A,F; 304 IS is_A_all; 305 PetscInt *idxs_schur,n_I; 306 307 /* get working mat */ 308 ierr = ISGetLocalSize(sub_schurs->is_I_layer,&n_I);CHKERRQ(ierr); 309 ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size+n_I,all_local_idx_N,PETSC_COPY_VALUES,&is_A_all);CHKERRQ(ierr); 310 ierr = MatGetSubMatrixUnsorted(sub_schurs->A,is_A_all,is_A_all,MAT_INITIAL_MATRIX,&A);CHKERRQ(ierr); 311 ierr = ISDestroy(&is_A_all);CHKERRQ(ierr); 312 313 if (n_I) { 314 if (is_symmetric) { 315 ierr = MatGetFactor(A,MATSOLVERMUMPS,MAT_FACTOR_CHOLESKY,&F);CHKERRQ(ierr); 316 } else { 317 ierr = MatGetFactor(A,MATSOLVERMUMPS,MAT_FACTOR_LU,&F);CHKERRQ(ierr); 318 } 319 320 /* subsets ordered last */ 321 ierr = PetscMalloc1(local_size,&idxs_schur);CHKERRQ(ierr); 322 for (i=0;i<local_size;i++) { 323 idxs_schur[i] = n_I+i+1; 324 } 325 ierr = MatMumpsSetSchurIndices(F,local_size,idxs_schur);CHKERRQ(ierr); 326 ierr = PetscFree(idxs_schur);CHKERRQ(ierr); 327 328 /* factorization step */ 329 if (is_symmetric) { 330 ierr = MatCholeskyFactorSymbolic(F,A,NULL,NULL);CHKERRQ(ierr); 331 ierr = MatCholeskyFactorNumeric(F,A,NULL);CHKERRQ(ierr); 332 } else { 333 ierr = MatLUFactorSymbolic(F,A,NULL,NULL,NULL);CHKERRQ(ierr); 334 ierr = MatLUFactorNumeric(F,A,NULL);CHKERRQ(ierr); 335 } 336 337 /* get explicit Schur Complement computed during numeric factorization */ 338 ierr = MatMumpsGetSchurComplement(F,&S_all);CHKERRQ(ierr); 339 ierr = MatDestroy(&F);CHKERRQ(ierr); 340 } else { 341 ierr = MatConvert(A,MATSEQDENSE,MAT_INITIAL_MATRIX,&S_all);CHKERRQ(ierr); 342 } 343 ierr = MatDestroy(&A);CHKERRQ(ierr); 344 #endif 345 } 346 347 if (!sub_schurs->n_subs_seq_g) { 348 PetscFunctionReturn(0); 349 } 350 351 /* subset indices in local boundary numbering */ 352 if (!sub_schurs->is_Ej_all) { 353 ierr = ISGlobalToLocalMappingApply(sub_schurs->BtoNmap,IS_GTOLM_DROP,local_size,all_local_idx_N+extra,&subset_size,all_local_idx_B);CHKERRQ(ierr); 354 if (subset_size != local_size) { 355 SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in sub_schurs serial (BtoNmap)! %d != %d\n",subset_size,local_size); 356 } 357 ierr = ISCreateGeneral(PETSC_COMM_SELF,local_size,all_local_idx_B,PETSC_OWN_POINTER,&sub_schurs->is_Ej_all);CHKERRQ(ierr); 358 } 359 360 /* Local matrix of all local Schur on subsets transposed */ 361 if (!sub_schurs->S_Ej_all) { 362 ierr = MatCreate(PETSC_COMM_SELF,&sub_schurs->S_Ej_all);CHKERRQ(ierr); 363 ierr = MatSetSizes(sub_schurs->S_Ej_all,PETSC_DECIDE,PETSC_DECIDE,local_size,local_size);CHKERRQ(ierr); 364 ierr = MatSetType(sub_schurs->S_Ej_all,MATAIJ);CHKERRQ(ierr); 365 ierr = MatSeqAIJSetPreallocation(sub_schurs->S_Ej_all,0,nnz);CHKERRQ(ierr); 366 } else { 367 ierr = MatZeroEntries(sub_schurs->S_Ej_all);CHKERRQ(ierr); 368 } 369 370 S_Ej_tilda_all = 0; 371 S_Ej_inv_all = 0; 372 work = NULL; 373 pivots = NULL; 374 Stilda_computed = PETSC_FALSE; 375 if (sub_schurs->n_subs && deluxe) { /* workspace needed only for GETRI */ 376 PetscScalar lwork; 377 378 B_lwork = -1; 379 ierr = PetscBLASIntCast(max_subset_size,&B_N);CHKERRQ(ierr); 380 ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 381 PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,work,&B_N,pivots,&lwork,&B_lwork,&B_ierr)); 382 ierr = PetscFPTrapPop();CHKERRQ(ierr); 383 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in query to GETRI Lapack routine %d",(int)B_ierr); 384 ierr = PetscBLASIntCast((PetscInt)lwork,&B_lwork);CHKERRQ(ierr); 385 ierr = PetscMalloc2(B_lwork,&work,max_subset_size,&pivots);CHKERRQ(ierr); 386 } 387 388 ierr = PetscBTMemzero(sub_schurs->n_subs,sub_schurs->computed_Stilda_subs);CHKERRQ(ierr); 389 if (!sub_schurs->use_mumps) { 390 PetscScalar *fill_vals; 391 PetscInt *dummy_idx; 392 393 /* Work arrays */ 394 ierr = PetscMalloc1(max_subset_size,&dummy_idx);CHKERRQ(ierr); 395 /* Loop on local problems end compute Schur complements explicitly */ 396 local_size = 0; 397 for (i=0;i<sub_schurs->n_subs_seq;i++) { 398 Mat S_Ej_expl; 399 PetscInt j,lpi = sub_schurs->index_sequential[i]; 400 PetscBool Sdense; 401 402 ierr = PCBDDCComputeExplicitSchur(sub_schurs->S_Ej[lpi],&S_Ej_expl);CHKERRQ(ierr); 403 ierr = ISGetLocalSize(sub_schurs->is_subs[lpi],&subset_size);CHKERRQ(ierr); 404 ierr = PetscObjectTypeCompare((PetscObject)S_Ej_expl,MATSEQDENSE,&Sdense);CHKERRQ(ierr); 405 if (Sdense) { 406 for (j=0;j<subset_size;j++) { 407 dummy_idx[j]=local_size+j; 408 } 409 ierr = MatDenseGetArray(S_Ej_expl,&fill_vals);CHKERRQ(ierr); 410 ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,fill_vals,INSERT_VALUES);CHKERRQ(ierr); 411 ierr = MatDenseRestoreArray(S_Ej_expl,&fill_vals);CHKERRQ(ierr); 412 } else { 413 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not yet implemented for sparse matrices"); 414 } 415 local_size += subset_size; 416 ierr = MatDestroy(&sub_schurs->S_Ej[lpi]);CHKERRQ(ierr); 417 sub_schurs->S_Ej[lpi] = S_Ej_expl; 418 } 419 /* Stildas are not computed without mumps */ 420 for (i=0;i<sub_schurs->n_subs;i++) { 421 ierr = MatDestroy(&sub_schurs->S_Ej_tilda[i]);CHKERRQ(ierr); 422 } 423 ierr = PetscFree(dummy_idx);CHKERRQ(ierr); 424 } else { 425 PetscInt *dummy_idx,n_all; 426 427 if (compute_Stilda) { 428 ierr = MatCreate(PETSC_COMM_SELF,&S_Ej_tilda_all);CHKERRQ(ierr); 429 ierr = MatSetSizes(S_Ej_tilda_all,PETSC_DECIDE,PETSC_DECIDE,local_size,local_size);CHKERRQ(ierr); 430 ierr = MatSetType(S_Ej_tilda_all,MATAIJ);CHKERRQ(ierr); 431 ierr = MatSeqAIJSetPreallocation(S_Ej_tilda_all,0,nnz);CHKERRQ(ierr); 432 if (deluxe) { 433 ierr = MatCreate(PETSC_COMM_SELF,&S_Ej_inv_all);CHKERRQ(ierr); 434 ierr = MatSetSizes(S_Ej_inv_all,PETSC_DECIDE,PETSC_DECIDE,local_size,local_size);CHKERRQ(ierr); 435 ierr = MatSetType(S_Ej_inv_all,MATAIJ);CHKERRQ(ierr); 436 ierr = MatSeqAIJSetPreallocation(S_Ej_inv_all,0,nnz);CHKERRQ(ierr); 437 } 438 } 439 440 ierr = MatGetSize(S_all,&n_all,NULL);CHKERRQ(ierr); 441 local_size = 0; 442 /* Work arrays */ 443 ierr = PetscMalloc1(max_subset_size,&dummy_idx);CHKERRQ(ierr); 444 445 Stilda_computed = PETSC_FALSE; 446 for (i=0;i<sub_schurs->n_subs;i++) { 447 IS is_E; 448 PetscScalar *vals; 449 PetscInt j; 450 451 ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr); 452 for (j=0;j<subset_size;j++) { 453 dummy_idx[j]=local_size+j; 454 } 455 ierr = ISCreateStride(PetscObjectComm((PetscObject)sub_schurs->is_I),subset_size,local_size,1,&is_E);CHKERRQ(ierr); 456 if (sub_schurs->S_Ej[i]) { 457 ierr = MatGetSubMatrix(S_all,is_E,is_E,MAT_REUSE_MATRIX,&sub_schurs->S_Ej[i]);CHKERRQ(ierr); 458 } else { 459 ierr = MatGetSubMatrix(S_all,is_E,is_E,MAT_INITIAL_MATRIX,&sub_schurs->S_Ej[i]);CHKERRQ(ierr); 460 } 461 ierr = MatDenseGetArray(sub_schurs->S_Ej[i],&vals);CHKERRQ(ierr); 462 ierr = MatSetValues(sub_schurs->S_Ej_all,subset_size,dummy_idx,subset_size,dummy_idx,vals,INSERT_VALUES);CHKERRQ(ierr); 463 ierr = MatDenseRestoreArray(sub_schurs->S_Ej[i],&vals);CHKERRQ(ierr); 464 465 if (compute_Stilda && ((PetscBTLookup(sub_schurs->is_edge,i) && use_edges) || use_faces)) { 466 ierr = PetscBTSet(sub_schurs->computed_Stilda_subs,i);CHKERRQ(ierr); 467 Stilda_computed = PETSC_TRUE; 468 } 469 470 if (PetscBTLookup(sub_schurs->computed_Stilda_subs,i)) { 471 Mat Stilda; 472 if (sub_schurs->n_subs == 1) { 473 ierr = PetscObjectReference((PetscObject)sub_schurs->S_Ej[i]);CHKERRQ(ierr); 474 Stilda = sub_schurs->S_Ej[i]; 475 } else { 476 KSP ksp; 477 PC pc; 478 Mat S_EF,S_FE,S_FF,Stilda_impl; 479 IS is_F; 480 PetscScalar eps=1.e-8; 481 PetscBool chop=PETSC_FALSE; 482 483 ierr = ISComplement(is_E,0,n_all,&is_F);CHKERRQ(ierr); 484 ierr = MatGetSubMatrix(S_all,is_E,is_F,MAT_INITIAL_MATRIX,&S_EF);CHKERRQ(ierr); 485 ierr = MatGetSubMatrix(S_all,is_F,is_F,MAT_INITIAL_MATRIX,&S_FF);CHKERRQ(ierr); 486 ierr = MatGetSubMatrix(S_all,is_F,is_E,MAT_INITIAL_MATRIX,&S_FE);CHKERRQ(ierr); 487 ierr = ISDestroy(&is_F);CHKERRQ(ierr); 488 if (chop) { 489 ierr = MatChop(S_FF,eps);CHKERRQ(ierr); 490 ierr = MatConvert(S_FF,MATAIJ,MAT_REUSE_MATRIX,&S_FF);CHKERRQ(ierr); 491 } 492 ierr = MatCreateSchurComplement(S_FF,S_FF,S_FE,S_EF,sub_schurs->S_Ej[i],&Stilda_impl);CHKERRQ(ierr); 493 ierr = MatDestroy(&S_FF);CHKERRQ(ierr); 494 ierr = MatDestroy(&S_FE);CHKERRQ(ierr); 495 ierr = MatDestroy(&S_EF);CHKERRQ(ierr); 496 ierr = MatSchurComplementGetKSP(Stilda_impl,&ksp);CHKERRQ(ierr); 497 ierr = KSPSetType(ksp,KSPPREONLY);CHKERRQ(ierr); 498 ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr); 499 if (is_symmetric) { 500 ierr = PCSetType(pc,PCCHOLESKY);CHKERRQ(ierr); 501 } else { 502 ierr = PCSetType(pc,PCLU);CHKERRQ(ierr); 503 } 504 if (!chop) { 505 ierr = PCFactorSetUseInPlace(pc,PETSC_TRUE);CHKERRQ(ierr); 506 } else { 507 ierr = PCFactorSetMatSolverPackage(pc,MATSOLVERMUMPS);CHKERRQ(ierr); 508 } 509 ierr = KSPSetUp(ksp);CHKERRQ(ierr); 510 ierr = PCBDDCComputeExplicitSchur(Stilda_impl,&Stilda);CHKERRQ(ierr); 511 ierr = MatDestroy(&Stilda_impl);CHKERRQ(ierr); 512 } 513 /* 514 PetscViewerPushFormat(PETSC_VIEWER_STDOUT_SELF,PETSC_VIEWER_ASCII_MATLAB); 515 ierr = MatView(Stilda,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); 516 PetscViewerPopFormat(PETSC_VIEWER_STDOUT_SELF); 517 */ 518 ierr = MatDestroy(&sub_schurs->S_Ej_tilda[i]);CHKERRQ(ierr); 519 sub_schurs->S_Ej_tilda[i] = Stilda; 520 521 522 ierr = MatDenseGetArray(sub_schurs->S_Ej_tilda[i],&vals);CHKERRQ(ierr); 523 if (deluxe) { /* when using deluxe scaling, we need (S_1^-1+S_2^-1)^-1 */ 524 PetscScalar *vals2; 525 526 ierr = MatDenseGetArray(sub_schurs->S_Ej[i],&vals2);CHKERRQ(ierr); 527 /* need to be optimized (cholesky) */ 528 ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr); 529 ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 530 PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,vals,&B_N,pivots,&B_ierr)); 531 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr); 532 PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,vals,&B_N,pivots,work,&B_lwork,&B_ierr)); 533 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr); 534 PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,vals2,&B_N,pivots,&B_ierr)); 535 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr); 536 PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,vals2,&B_N,pivots,work,&B_lwork,&B_ierr)); 537 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr); 538 ierr = PetscFPTrapPop();CHKERRQ(ierr); 539 ierr = MatSetValues(S_Ej_inv_all,subset_size,dummy_idx,subset_size,dummy_idx,vals2,INSERT_VALUES);CHKERRQ(ierr); 540 ierr = MatDenseRestoreArray(sub_schurs->S_Ej[i],&vals2);CHKERRQ(ierr); 541 } 542 ierr = MatSetValues(S_Ej_tilda_all,subset_size,dummy_idx,subset_size,dummy_idx,vals,INSERT_VALUES);CHKERRQ(ierr); 543 ierr = MatDenseRestoreArray(sub_schurs->S_Ej_tilda[i],&vals);CHKERRQ(ierr); 544 } else { 545 ierr = MatDestroy(&sub_schurs->S_Ej_tilda[i]);CHKERRQ(ierr); 546 } 547 ierr = ISDestroy(&is_E);CHKERRQ(ierr); 548 local_size += subset_size; 549 } 550 ierr = PetscFree(dummy_idx);CHKERRQ(ierr); 551 } 552 ierr = PetscFree(nnz);CHKERRQ(ierr); 553 ierr = MatDestroy(&S_all);CHKERRQ(ierr); 554 ierr = MatAssemblyBegin(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 555 ierr = MatAssemblyEnd(sub_schurs->S_Ej_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 556 if (compute_Stilda) { 557 ierr = MatAssemblyBegin(S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 558 ierr = MatAssemblyEnd(S_Ej_tilda_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 559 if (deluxe) { 560 ierr = MatAssemblyBegin(S_Ej_inv_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 561 ierr = MatAssemblyEnd(S_Ej_inv_all,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 562 } 563 } 564 565 /* Global matrix of all assembled Schur on subsets */ 566 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); 567 ierr = ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)sub_schurs->l2gmap),1,local_size,all_local_idx_G,PETSC_COPY_VALUES,&l2gmap_subsets);CHKERRQ(ierr); 568 ierr = MatCreateIS(PetscObjectComm((PetscObject)sub_schurs->l2gmap),1,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size,l2gmap_subsets,&work_mat);CHKERRQ(ierr); 569 ierr = ISLocalToGlobalMappingDestroy(&l2gmap_subsets);CHKERRQ(ierr); 570 ierr = MatISSetLocalMat(work_mat,sub_schurs->S_Ej_all);CHKERRQ(ierr); 571 ierr = MatCreate(PetscObjectComm((PetscObject)work_mat),&global_schur_subsets);CHKERRQ(ierr); 572 ierr = MatSetSizes(global_schur_subsets,PETSC_DECIDE,PETSC_DECIDE,global_size,global_size);CHKERRQ(ierr); 573 ierr = MatSetType(global_schur_subsets,MATMPIAIJ);CHKERRQ(ierr); 574 ierr = MatISSetMPIXAIJPreallocation_Private(work_mat,global_schur_subsets,PETSC_TRUE);CHKERRQ(ierr); 575 ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr); 576 /* Get local part of (\sum_j S_Ej) */ 577 ierr = PetscFree(all_local_idx_N);CHKERRQ(ierr); 578 ierr = ISCreateGeneral(PetscObjectComm((PetscObject)sub_schurs->l2gmap),local_size,all_local_idx_G,PETSC_OWN_POINTER,&temp_is);CHKERRQ(ierr); 579 if (sub_schurs->sum_S_Ej_all) { 580 ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,MAT_REUSE_MATRIX,&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr); 581 } else { 582 ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,MAT_INITIAL_MATRIX,&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr); 583 } 584 585 if (Stilda_computed) { 586 ierr = MatISSetLocalMat(work_mat,S_Ej_tilda_all);CHKERRQ(ierr); 587 ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr); 588 if (sub_schurs->sum_S_Ej_tilda_all) { 589 ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,MAT_REUSE_MATRIX,&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr); 590 } else { 591 ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,MAT_INITIAL_MATRIX,&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr); 592 } 593 if (deluxe) { 594 PetscInt cum; 595 PetscScalar *array,*array2; 596 ierr = MatISSetLocalMat(work_mat,S_Ej_inv_all);CHKERRQ(ierr); 597 ierr = MatISGetMPIXAIJ(work_mat,MAT_REUSE_MATRIX,&global_schur_subsets);CHKERRQ(ierr); 598 if (sub_schurs->sum_S_Ej_inv_all) { 599 ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,MAT_REUSE_MATRIX,&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr); 600 } else { 601 ierr = MatGetSubMatrixUnsorted(global_schur_subsets,temp_is,temp_is,MAT_INITIAL_MATRIX,&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr); 602 } 603 /* invert blocks of sum_S_Ej_inv_all */ 604 ierr = MatSeqAIJGetArray(sub_schurs->sum_S_Ej_inv_all,&array);CHKERRQ(ierr); 605 ierr = MatSeqAIJGetArray(sub_schurs->sum_S_Ej_tilda_all,&array2);CHKERRQ(ierr); 606 cum = 0; 607 for (i=0;i<sub_schurs->n_subs;i++) { 608 ierr = ISGetLocalSize(sub_schurs->is_subs[i],&subset_size);CHKERRQ(ierr); 609 if (PetscBTLookup(sub_schurs->computed_Stilda_subs,i)) { 610 /* need to be optimized (cholesky) */ 611 ierr = PetscBLASIntCast(subset_size,&B_N);CHKERRQ(ierr); 612 ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 613 PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,array+cum,&B_N,pivots,&B_ierr)); 614 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr); 615 PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,array+cum,&B_N,pivots,work,&B_lwork,&B_ierr)); 616 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr); 617 PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,array2+cum,&B_N,pivots,&B_ierr)); 618 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr); 619 PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,array2+cum,&B_N,pivots,work,&B_lwork,&B_ierr)); 620 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr); 621 ierr = PetscFPTrapPop();CHKERRQ(ierr); 622 } 623 cum += subset_size*subset_size; 624 } 625 ierr = MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_inv_all,&array);CHKERRQ(ierr); 626 ierr = MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_tilda_all,&array2);CHKERRQ(ierr); 627 } 628 } 629 630 ierr = PetscFree2(work,pivots);CHKERRQ(ierr); 631 ierr = MatDestroy(&global_schur_subsets);CHKERRQ(ierr); 632 ierr = MatDestroy(&S_Ej_tilda_all);CHKERRQ(ierr); 633 ierr = MatDestroy(&S_Ej_inv_all);CHKERRQ(ierr); 634 ierr = MatDestroy(&work_mat);CHKERRQ(ierr); 635 ierr = ISDestroy(&temp_is);CHKERRQ(ierr); 636 PetscFunctionReturn(0); 637 } 638 639 #undef __FUNCT__ 640 #define __FUNCT__ "PCBDDCSubSchursInit" 641 PetscErrorCode PCBDDCSubSchursInit(PCBDDCSubSchurs sub_schurs, Mat A, Mat S, IS is_I, IS is_B, PCBDDCGraph graph, ISLocalToGlobalMapping BtoNmap, PetscInt seqthreshold) 642 { 643 IS *faces,*edges,*all_cc,vertices; 644 PetscInt *index_sequential,*index_parallel; 645 PetscInt *auxlocal_sequential,*auxlocal_parallel; 646 PetscInt *auxglobal_sequential,*auxglobal_parallel; 647 PetscInt *auxmapping; 648 PetscInt i,max_subset_size; 649 PetscInt n_sequential_problems,n_local_sequential_problems,n_parallel_problems,n_local_parallel_problems; 650 PetscInt n_faces,n_edges,n_all_cc; 651 PetscBool is_sorted; 652 PetscErrorCode ierr; 653 654 PetscFunctionBegin; 655 ierr = ISSorted(is_I,&is_sorted);CHKERRQ(ierr); 656 if (!is_sorted) { 657 SETERRQ(PetscObjectComm((PetscObject)is_I),PETSC_ERR_PLIB,"IS for I dofs should be shorted"); 658 } 659 ierr = ISSorted(is_B,&is_sorted);CHKERRQ(ierr); 660 if (!is_sorted) { 661 SETERRQ(PetscObjectComm((PetscObject)is_B),PETSC_ERR_PLIB,"IS for B dofs should be shorted"); 662 } 663 664 /* reset any previous data */ 665 ierr = PCBDDCSubSchursReset(sub_schurs);CHKERRQ(ierr); 666 667 /* get index sets for faces and edges */ 668 ierr = PCBDDCGraphGetCandidatesIS(graph,&n_faces,&faces,&n_edges,&edges,&vertices);CHKERRQ(ierr); 669 n_all_cc = n_faces+n_edges; 670 ierr = PetscBTCreate(n_all_cc,&sub_schurs->is_edge);CHKERRQ(ierr); 671 ierr = PetscBTCreate(n_all_cc,&sub_schurs->computed_Stilda_subs);CHKERRQ(ierr); 672 ierr = PetscMalloc1(n_all_cc,&all_cc);CHKERRQ(ierr); 673 for (i=0;i<n_faces;i++) { 674 all_cc[i] = faces[i]; 675 } 676 for (i=0;i<n_edges;i++) { 677 all_cc[n_faces+i] = edges[i]; 678 ierr = PetscBTSet(sub_schurs->is_edge,n_faces+i);CHKERRQ(ierr); 679 } 680 ierr = PetscFree(faces);CHKERRQ(ierr); 681 ierr = PetscFree(edges);CHKERRQ(ierr); 682 683 /* map interface's subsets */ 684 max_subset_size = 0; 685 for (i=0;i<n_all_cc;i++) { 686 PetscInt subset_size; 687 ierr = ISGetLocalSize(all_cc[i],&subset_size);CHKERRQ(ierr); 688 max_subset_size = PetscMax(max_subset_size,subset_size); 689 } 690 ierr = PetscMalloc1(max_subset_size,&auxmapping);CHKERRQ(ierr); 691 ierr = PetscMalloc2(graph->ncc,&auxlocal_sequential, 692 graph->ncc,&auxlocal_parallel);CHKERRQ(ierr); 693 ierr = PetscMalloc2(graph->ncc,&index_sequential, 694 graph->ncc,&index_parallel);CHKERRQ(ierr); 695 696 /* if threshold is negative uses all sequential problems (possibly using MUMPS) */ 697 sub_schurs->use_mumps = PETSC_FALSE; 698 if (seqthreshold < 0) { 699 seqthreshold = max_subset_size; 700 #if defined(PETSC_HAVE_MUMPS) 701 sub_schurs->use_mumps = !!A; 702 #endif 703 } 704 705 706 /* determine which problem has to be solved in parallel or sequentially */ 707 n_local_sequential_problems = 0; 708 n_local_parallel_problems = 0; 709 for (i=0;i<n_all_cc;i++) { 710 PetscInt subset_size,j,min_loc = 0; 711 const PetscInt *idxs; 712 713 ierr = ISGetLocalSize(all_cc[i],&subset_size);CHKERRQ(ierr); 714 ierr = ISGetIndices(all_cc[i],&idxs);CHKERRQ(ierr); 715 ierr = ISLocalToGlobalMappingApply(graph->l2gmap,subset_size,idxs,auxmapping);CHKERRQ(ierr); 716 for (j=1;j<subset_size;j++) { 717 if (auxmapping[j]<auxmapping[min_loc]) { 718 min_loc = j; 719 } 720 } 721 if (subset_size > seqthreshold) { 722 index_parallel[n_local_parallel_problems] = i; 723 auxlocal_parallel[n_local_parallel_problems] = idxs[min_loc]; 724 n_local_parallel_problems++; 725 } else { 726 index_sequential[n_local_sequential_problems] = i; 727 auxlocal_sequential[n_local_sequential_problems] = idxs[min_loc]; 728 n_local_sequential_problems++; 729 } 730 ierr = ISRestoreIndices(all_cc[i],&idxs);CHKERRQ(ierr); 731 } 732 733 /* Number parallel problems */ 734 auxglobal_parallel = 0; 735 ierr = PCBDDCSubsetNumbering(PetscObjectComm((PetscObject)graph->l2gmap),graph->l2gmap,n_local_parallel_problems,auxlocal_parallel,PETSC_NULL,&n_parallel_problems,&auxglobal_parallel);CHKERRQ(ierr); 736 737 /* Number sequential problems */ 738 auxglobal_sequential = 0; 739 ierr = PCBDDCSubsetNumbering(PetscObjectComm((PetscObject)graph->l2gmap),graph->l2gmap,n_local_sequential_problems,auxlocal_sequential,PETSC_NULL,&n_sequential_problems,&auxglobal_sequential);CHKERRQ(ierr); 740 741 /* update info in sub_schurs */ 742 if (sub_schurs->use_mumps && A) { 743 ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 744 sub_schurs->A = A; 745 } 746 ierr = PetscObjectReference((PetscObject)S);CHKERRQ(ierr); 747 sub_schurs->S = S; 748 ierr = PetscObjectReference((PetscObject)is_I);CHKERRQ(ierr); 749 sub_schurs->is_I = is_I; 750 ierr = PetscObjectReference((PetscObject)is_B);CHKERRQ(ierr); 751 sub_schurs->is_B = is_B; 752 ierr = PetscObjectReference((PetscObject)graph->l2gmap);CHKERRQ(ierr); 753 sub_schurs->l2gmap = graph->l2gmap; 754 ierr = PetscObjectReference((PetscObject)BtoNmap);CHKERRQ(ierr); 755 sub_schurs->BtoNmap = BtoNmap; 756 sub_schurs->n_subs_seq = n_local_sequential_problems; 757 sub_schurs->n_subs_par = n_local_parallel_problems; 758 sub_schurs->n_subs_seq_g = n_sequential_problems; 759 sub_schurs->n_subs_par_g = n_parallel_problems; 760 sub_schurs->n_subs = sub_schurs->n_subs_seq + sub_schurs->n_subs_par; 761 sub_schurs->is_subs = all_cc; 762 if (!sub_schurs->use_mumps) { /* for adaptive selection */ 763 for (i=0;i<sub_schurs->n_subs;i++) { 764 ierr = ISSort(sub_schurs->is_subs[i]);CHKERRQ(ierr); 765 } 766 } 767 sub_schurs->is_Ej_com = vertices; 768 sub_schurs->index_sequential = index_sequential; 769 sub_schurs->index_parallel = index_parallel; 770 sub_schurs->auxglobal_sequential = auxglobal_sequential; 771 sub_schurs->auxglobal_parallel = auxglobal_parallel; 772 773 774 /* allocate space for schur complements */ 775 ierr = PetscCalloc2(sub_schurs->n_subs,&sub_schurs->S_Ej, 776 sub_schurs->n_subs,&sub_schurs->S_Ej_tilda);CHKERRQ(ierr); 777 sub_schurs->S_Ej_all = NULL; 778 sub_schurs->sum_S_Ej_all = NULL; 779 sub_schurs->sum_S_Ej_inv_all = NULL; 780 sub_schurs->sum_S_Ej_tilda_all = NULL; 781 sub_schurs->is_Ej_all = NULL; 782 783 /* free workspace */ 784 ierr = PetscFree(auxmapping);CHKERRQ(ierr); 785 ierr = PetscFree2(auxlocal_sequential,auxlocal_parallel);CHKERRQ(ierr); 786 787 PetscFunctionReturn(0); 788 } 789 790 #undef __FUNCT__ 791 #define __FUNCT__ "PCBDDCSubSchursCreate" 792 PetscErrorCode PCBDDCSubSchursCreate(PCBDDCSubSchurs *sub_schurs) 793 { 794 PCBDDCSubSchurs schurs_ctx; 795 PetscErrorCode ierr; 796 797 PetscFunctionBegin; 798 ierr = PetscNew(&schurs_ctx);CHKERRQ(ierr); 799 schurs_ctx->n_subs = 0; 800 *sub_schurs = schurs_ctx; 801 PetscFunctionReturn(0); 802 } 803 804 #undef __FUNCT__ 805 #define __FUNCT__ "PCBDDCSubSchursDestroy" 806 PetscErrorCode PCBDDCSubSchursDestroy(PCBDDCSubSchurs *sub_schurs) 807 { 808 PetscErrorCode ierr; 809 810 PetscFunctionBegin; 811 ierr = PCBDDCSubSchursReset(*sub_schurs);CHKERRQ(ierr); 812 ierr = PetscFree(*sub_schurs);CHKERRQ(ierr); 813 PetscFunctionReturn(0); 814 } 815 816 #undef __FUNCT__ 817 #define __FUNCT__ "PCBDDCSubSchursReset" 818 PetscErrorCode PCBDDCSubSchursReset(PCBDDCSubSchurs sub_schurs) 819 { 820 PetscInt i; 821 PetscErrorCode ierr; 822 823 PetscFunctionBegin; 824 ierr = MatDestroy(&sub_schurs->A);CHKERRQ(ierr); 825 ierr = MatDestroy(&sub_schurs->S);CHKERRQ(ierr); 826 ierr = ISDestroy(&sub_schurs->is_I);CHKERRQ(ierr); 827 ierr = ISDestroy(&sub_schurs->is_B);CHKERRQ(ierr); 828 ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->l2gmap);CHKERRQ(ierr); 829 ierr = ISLocalToGlobalMappingDestroy(&sub_schurs->BtoNmap);CHKERRQ(ierr); 830 ierr = MatDestroy(&sub_schurs->S_Ej_all);CHKERRQ(ierr); 831 ierr = MatDestroy(&sub_schurs->sum_S_Ej_all);CHKERRQ(ierr); 832 ierr = MatDestroy(&sub_schurs->sum_S_Ej_inv_all);CHKERRQ(ierr); 833 ierr = MatDestroy(&sub_schurs->sum_S_Ej_tilda_all);CHKERRQ(ierr); 834 ierr = ISDestroy(&sub_schurs->is_Ej_all);CHKERRQ(ierr); 835 ierr = ISDestroy(&sub_schurs->is_Ej_com);CHKERRQ(ierr); 836 ierr = PetscBTDestroy(&sub_schurs->is_edge);CHKERRQ(ierr); 837 ierr = PetscBTDestroy(&sub_schurs->computed_Stilda_subs);CHKERRQ(ierr); 838 for (i=0;i<sub_schurs->n_subs;i++) { 839 ierr = ISDestroy(&sub_schurs->is_subs[i]);CHKERRQ(ierr); 840 ierr = MatDestroy(&sub_schurs->S_Ej[i]);CHKERRQ(ierr); 841 ierr = MatDestroy(&sub_schurs->S_Ej_tilda[i]);CHKERRQ(ierr); 842 } 843 ierr = ISDestroy(&sub_schurs->is_I_layer);CHKERRQ(ierr); 844 if (sub_schurs->n_subs) { 845 ierr = PetscFree(sub_schurs->is_subs);CHKERRQ(ierr); 846 ierr = PetscFree2(sub_schurs->S_Ej,sub_schurs->S_Ej_tilda);CHKERRQ(ierr); 847 ierr = PetscFree2(sub_schurs->index_sequential,sub_schurs->index_parallel);CHKERRQ(ierr); 848 ierr = PetscFree(sub_schurs->auxglobal_sequential);CHKERRQ(ierr); 849 ierr = PetscFree(sub_schurs->auxglobal_parallel);CHKERRQ(ierr); 850 } 851 sub_schurs->n_subs = 0; 852 PetscFunctionReturn(0); 853 } 854 855 #undef __FUNCT__ 856 #define __FUNCT__ "PCBDDCAdjGetNextLayer_Private" 857 PETSC_STATIC_INLINE PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt* queue_tip,PetscInt n_prev,PetscBT touched,PetscInt* xadj,PetscInt* adjncy,PetscInt* n_added) 858 { 859 PetscInt i,j,n; 860 PetscErrorCode ierr; 861 862 PetscFunctionBegin; 863 n = 0; 864 for (i=-n_prev;i<0;i++) { 865 PetscInt start_dof = queue_tip[i]; 866 for (j=xadj[start_dof];j<xadj[start_dof+1];j++) { 867 PetscInt dof = adjncy[j]; 868 if (!PetscBTLookup(touched,dof)) { 869 ierr = PetscBTSet(touched,dof);CHKERRQ(ierr); 870 queue_tip[n] = dof; 871 n++; 872 } 873 } 874 } 875 *n_added = n; 876 PetscFunctionReturn(0); 877 } 878