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