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