1 #include <petsc/private/pcbddcimpl.h> 2 #include <petsc/private/pcbddcprivateimpl.h> 3 #include <../src/mat/impls/dense/seq/dense.h> 4 #include <petscblaslapack.h> 5 6 static inline PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt *, PetscInt, PetscBT, PetscInt *, PetscInt *, PetscInt *); 7 static PetscErrorCode PCBDDCComputeExplicitSchur(Mat, PetscBool, MatReuse, Mat *); 8 static PetscErrorCode PCBDDCReuseSolvers_Interior(PC, Vec, Vec); 9 static PetscErrorCode PCBDDCReuseSolvers_Correction(PC, Vec, Vec); 10 11 /* if v2 is not present, correction is done in-place */ 12 PetscErrorCode PCBDDCReuseSolversBenignAdapt(PCBDDCReuseSolvers ctx, Vec v, Vec v2, PetscBool sol, PetscBool full) 13 { 14 PetscScalar *array; 15 PetscScalar *array2; 16 17 PetscFunctionBegin; 18 if (!ctx->benign_n) PetscFunctionReturn(PETSC_SUCCESS); 19 if (sol && full) { 20 PetscInt n_I, size_schur; 21 22 /* get sizes */ 23 PetscCall(MatGetSize(ctx->benign_csAIB, &size_schur, NULL)); 24 PetscCall(VecGetSize(v, &n_I)); 25 n_I = n_I - size_schur; 26 /* get schur sol from array */ 27 PetscCall(VecGetArray(v, &array)); 28 PetscCall(VecPlaceArray(ctx->benign_dummy_schur_vec, array + n_I)); 29 PetscCall(VecRestoreArray(v, &array)); 30 /* apply interior sol correction */ 31 PetscCall(MatMultTranspose(ctx->benign_csAIB, ctx->benign_dummy_schur_vec, ctx->benign_corr_work)); 32 PetscCall(VecResetArray(ctx->benign_dummy_schur_vec)); 33 PetscCall(MatMultAdd(ctx->benign_AIIm1ones, ctx->benign_corr_work, v, v)); 34 } 35 if (v2) { 36 PetscInt nl; 37 38 PetscCall(VecGetArrayRead(v, (const PetscScalar **)&array)); 39 PetscCall(VecGetLocalSize(v2, &nl)); 40 PetscCall(VecGetArray(v2, &array2)); 41 PetscCall(PetscArraycpy(array2, array, nl)); 42 } else { 43 PetscCall(VecGetArray(v, &array)); 44 array2 = array; 45 } 46 if (!sol) { /* change rhs */ 47 PetscInt n; 48 for (n = 0; n < ctx->benign_n; n++) { 49 PetscScalar sum = 0.; 50 const PetscInt *cols; 51 PetscInt nz, i; 52 53 PetscCall(ISGetLocalSize(ctx->benign_zerodiag_subs[n], &nz)); 54 PetscCall(ISGetIndices(ctx->benign_zerodiag_subs[n], &cols)); 55 for (i = 0; i < nz - 1; i++) sum += array[cols[i]]; 56 #if defined(PETSC_USE_COMPLEX) 57 sum = -(PetscRealPart(sum) / nz + PETSC_i * (PetscImaginaryPart(sum) / nz)); 58 #else 59 sum = -sum / nz; 60 #endif 61 for (i = 0; i < nz - 1; i++) array2[cols[i]] += sum; 62 ctx->benign_save_vals[n] = array2[cols[nz - 1]]; 63 array2[cols[nz - 1]] = sum; 64 PetscCall(ISRestoreIndices(ctx->benign_zerodiag_subs[n], &cols)); 65 } 66 } else { 67 PetscInt n; 68 for (n = 0; n < ctx->benign_n; n++) { 69 PetscScalar sum = 0.; 70 const PetscInt *cols; 71 PetscInt nz, i; 72 PetscCall(ISGetLocalSize(ctx->benign_zerodiag_subs[n], &nz)); 73 PetscCall(ISGetIndices(ctx->benign_zerodiag_subs[n], &cols)); 74 for (i = 0; i < nz - 1; i++) sum += array[cols[i]]; 75 #if defined(PETSC_USE_COMPLEX) 76 sum = -(PetscRealPart(sum) / nz + PETSC_i * (PetscImaginaryPart(sum) / nz)); 77 #else 78 sum = -sum / nz; 79 #endif 80 for (i = 0; i < nz - 1; i++) array2[cols[i]] += sum; 81 array2[cols[nz - 1]] = ctx->benign_save_vals[n]; 82 PetscCall(ISRestoreIndices(ctx->benign_zerodiag_subs[n], &cols)); 83 } 84 } 85 if (v2) { 86 PetscCall(VecRestoreArrayRead(v, (const PetscScalar **)&array)); 87 PetscCall(VecRestoreArray(v2, &array2)); 88 } else { 89 PetscCall(VecRestoreArray(v, &array)); 90 } 91 if (!sol && full) { 92 Vec usedv; 93 PetscInt n_I, size_schur; 94 95 /* get sizes */ 96 PetscCall(MatGetSize(ctx->benign_csAIB, &size_schur, NULL)); 97 PetscCall(VecGetSize(v, &n_I)); 98 n_I = n_I - size_schur; 99 /* compute schur rhs correction */ 100 if (v2) { 101 usedv = v2; 102 } else { 103 usedv = v; 104 } 105 /* apply schur rhs correction */ 106 PetscCall(MatMultTranspose(ctx->benign_AIIm1ones, usedv, ctx->benign_corr_work)); 107 PetscCall(VecGetArrayRead(usedv, (const PetscScalar **)&array)); 108 PetscCall(VecPlaceArray(ctx->benign_dummy_schur_vec, array + n_I)); 109 PetscCall(VecRestoreArrayRead(usedv, (const PetscScalar **)&array)); 110 PetscCall(MatMultAdd(ctx->benign_csAIB, ctx->benign_corr_work, ctx->benign_dummy_schur_vec, ctx->benign_dummy_schur_vec)); 111 PetscCall(VecResetArray(ctx->benign_dummy_schur_vec)); 112 } 113 PetscFunctionReturn(PETSC_SUCCESS); 114 } 115 116 static PetscErrorCode PCBDDCReuseSolvers_Solve_Private(PC pc, Vec rhs, Vec sol, PetscBool transpose, PetscBool full) 117 { 118 PCBDDCReuseSolvers ctx; 119 PetscBool copy = PETSC_FALSE; 120 121 PetscFunctionBegin; 122 PetscCall(PCShellGetContext(pc, &ctx)); 123 if (full) { 124 PetscCall(MatMumpsSetIcntl(ctx->F, 26, -1)); 125 #if defined(PETSC_HAVE_MKL_PARDISO) 126 PetscCall(MatMkl_PardisoSetCntl(ctx->F, 70, 0)); 127 #endif 128 copy = ctx->has_vertices; 129 } else { /* interior solver */ 130 PetscCall(MatMumpsSetIcntl(ctx->F, 26, 0)); 131 #if defined(PETSC_HAVE_MKL_PARDISO) 132 PetscCall(MatMkl_PardisoSetCntl(ctx->F, 70, 1)); 133 #endif 134 copy = PETSC_TRUE; 135 } 136 /* copy rhs into factored matrix workspace */ 137 if (copy) { 138 PetscInt n; 139 PetscScalar *array, *array_solver; 140 141 PetscCall(VecGetLocalSize(rhs, &n)); 142 PetscCall(VecGetArrayRead(rhs, (const PetscScalar **)&array)); 143 PetscCall(VecGetArray(ctx->rhs, &array_solver)); 144 PetscCall(PetscArraycpy(array_solver, array, n)); 145 PetscCall(VecRestoreArray(ctx->rhs, &array_solver)); 146 PetscCall(VecRestoreArrayRead(rhs, (const PetscScalar **)&array)); 147 148 PetscCall(PCBDDCReuseSolversBenignAdapt(ctx, ctx->rhs, NULL, PETSC_FALSE, full)); 149 if (transpose) { 150 PetscCall(MatSolveTranspose(ctx->F, ctx->rhs, ctx->sol)); 151 } else { 152 PetscCall(MatSolve(ctx->F, ctx->rhs, ctx->sol)); 153 } 154 PetscCall(PCBDDCReuseSolversBenignAdapt(ctx, ctx->sol, NULL, PETSC_TRUE, full)); 155 156 /* get back data to caller worskpace */ 157 PetscCall(VecGetArrayRead(ctx->sol, (const PetscScalar **)&array_solver)); 158 PetscCall(VecGetArray(sol, &array)); 159 PetscCall(PetscArraycpy(array, array_solver, n)); 160 PetscCall(VecRestoreArray(sol, &array)); 161 PetscCall(VecRestoreArrayRead(ctx->sol, (const PetscScalar **)&array_solver)); 162 } else { 163 if (ctx->benign_n) { 164 PetscCall(PCBDDCReuseSolversBenignAdapt(ctx, rhs, ctx->rhs, PETSC_FALSE, full)); 165 if (transpose) { 166 PetscCall(MatSolveTranspose(ctx->F, ctx->rhs, sol)); 167 } else { 168 PetscCall(MatSolve(ctx->F, ctx->rhs, sol)); 169 } 170 PetscCall(PCBDDCReuseSolversBenignAdapt(ctx, sol, NULL, PETSC_TRUE, full)); 171 } else { 172 if (transpose) { 173 PetscCall(MatSolveTranspose(ctx->F, rhs, sol)); 174 } else { 175 PetscCall(MatSolve(ctx->F, rhs, sol)); 176 } 177 } 178 } 179 /* restore defaults */ 180 PetscCall(MatMumpsSetIcntl(ctx->F, 26, -1)); 181 #if defined(PETSC_HAVE_MKL_PARDISO) 182 PetscCall(MatMkl_PardisoSetCntl(ctx->F, 70, 0)); 183 #endif 184 PetscFunctionReturn(PETSC_SUCCESS); 185 } 186 187 static PetscErrorCode PCBDDCReuseSolvers_Correction(PC pc, Vec rhs, Vec sol) 188 { 189 PetscFunctionBegin; 190 PetscCall(PCBDDCReuseSolvers_Solve_Private(pc, rhs, sol, PETSC_FALSE, PETSC_TRUE)); 191 PetscFunctionReturn(PETSC_SUCCESS); 192 } 193 194 static PetscErrorCode PCBDDCReuseSolvers_CorrectionTranspose(PC pc, Vec rhs, Vec sol) 195 { 196 PetscFunctionBegin; 197 PetscCall(PCBDDCReuseSolvers_Solve_Private(pc, rhs, sol, PETSC_TRUE, PETSC_TRUE)); 198 PetscFunctionReturn(PETSC_SUCCESS); 199 } 200 201 static PetscErrorCode PCBDDCReuseSolvers_Interior(PC pc, Vec rhs, Vec sol) 202 { 203 PetscFunctionBegin; 204 PetscCall(PCBDDCReuseSolvers_Solve_Private(pc, rhs, sol, PETSC_FALSE, PETSC_FALSE)); 205 PetscFunctionReturn(PETSC_SUCCESS); 206 } 207 208 static PetscErrorCode PCBDDCReuseSolvers_InteriorTranspose(PC pc, Vec rhs, Vec sol) 209 { 210 PetscFunctionBegin; 211 PetscCall(PCBDDCReuseSolvers_Solve_Private(pc, rhs, sol, PETSC_TRUE, PETSC_FALSE)); 212 PetscFunctionReturn(PETSC_SUCCESS); 213 } 214 215 static PetscErrorCode PCBDDCReuseSolvers_View(PC pc, PetscViewer viewer) 216 { 217 PCBDDCReuseSolvers ctx; 218 PetscBool isascii; 219 220 PetscFunctionBegin; 221 PetscCall(PCShellGetContext(pc, &ctx)); 222 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii)); 223 if (isascii) PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_INFO)); 224 PetscCall(MatView(ctx->F, viewer)); 225 if (isascii) PetscCall(PetscViewerPopFormat(viewer)); 226 PetscFunctionReturn(PETSC_SUCCESS); 227 } 228 229 static PetscErrorCode PCBDDCReuseSolversReset(PCBDDCReuseSolvers reuse) 230 { 231 PetscInt i; 232 233 PetscFunctionBegin; 234 PetscCall(MatDestroy(&reuse->F)); 235 PetscCall(VecDestroy(&reuse->sol)); 236 PetscCall(VecDestroy(&reuse->rhs)); 237 PetscCall(PCDestroy(&reuse->interior_solver)); 238 PetscCall(PCDestroy(&reuse->correction_solver)); 239 PetscCall(ISDestroy(&reuse->is_R)); 240 PetscCall(ISDestroy(&reuse->is_B)); 241 PetscCall(VecScatterDestroy(&reuse->correction_scatter_B)); 242 PetscCall(VecDestroy(&reuse->sol_B)); 243 PetscCall(VecDestroy(&reuse->rhs_B)); 244 for (i = 0; i < reuse->benign_n; i++) PetscCall(ISDestroy(&reuse->benign_zerodiag_subs[i])); 245 PetscCall(PetscFree(reuse->benign_zerodiag_subs)); 246 PetscCall(PetscFree(reuse->benign_save_vals)); 247 PetscCall(MatDestroy(&reuse->benign_csAIB)); 248 PetscCall(MatDestroy(&reuse->benign_AIIm1ones)); 249 PetscCall(VecDestroy(&reuse->benign_corr_work)); 250 PetscCall(VecDestroy(&reuse->benign_dummy_schur_vec)); 251 PetscFunctionReturn(PETSC_SUCCESS); 252 } 253 254 static PetscErrorCode PCBDDCReuseSolvers_Destroy(PC pc) 255 { 256 PCBDDCReuseSolvers ctx; 257 258 PetscFunctionBegin; 259 PetscCall(PCShellGetContext(pc, &ctx)); 260 PetscCall(PCBDDCReuseSolversReset(ctx)); 261 PetscCall(PetscFree(ctx)); 262 PetscCall(PCShellSetContext(pc, NULL)); 263 PetscFunctionReturn(PETSC_SUCCESS); 264 } 265 266 static PetscErrorCode PCBDDCComputeExplicitSchur(Mat M, PetscBool issym, MatReuse reuse, Mat *S) 267 { 268 Mat B, C, D, Bd, Cd, AinvBd; 269 KSP ksp; 270 PC pc; 271 PetscBool isLU, isILU, isCHOL, Bdense, Cdense; 272 PetscReal fill = 2.0; 273 PetscInt n_I; 274 PetscMPIInt size; 275 276 PetscFunctionBegin; 277 PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)M), &size)); 278 PetscCheck(size == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Not for parallel matrices"); 279 if (reuse == MAT_REUSE_MATRIX) { 280 PetscBool Sdense; 281 282 PetscCall(PetscObjectTypeCompare((PetscObject)*S, MATSEQDENSE, &Sdense)); 283 PetscCheck(Sdense, PetscObjectComm((PetscObject)M), PETSC_ERR_SUP, "S should dense"); 284 } 285 PetscCall(MatSchurComplementGetSubMatrices(M, NULL, NULL, &B, &C, &D)); 286 PetscCall(MatSchurComplementGetKSP(M, &ksp)); 287 PetscCall(KSPGetPC(ksp, &pc)); 288 PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCLU, &isLU)); 289 PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCILU, &isILU)); 290 PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCCHOLESKY, &isCHOL)); 291 PetscCall(PetscObjectTypeCompare((PetscObject)B, MATSEQDENSE, &Bdense)); 292 PetscCall(PetscObjectTypeCompare((PetscObject)C, MATSEQDENSE, &Cdense)); 293 PetscCall(MatGetSize(B, &n_I, NULL)); 294 if (n_I) { 295 if (!Bdense) { 296 PetscCall(MatConvert(B, MATSEQDENSE, MAT_INITIAL_MATRIX, &Bd)); 297 } else { 298 Bd = B; 299 } 300 301 if (isLU || isILU || isCHOL) { 302 Mat fact; 303 PetscCall(KSPSetUp(ksp)); 304 PetscCall(PCFactorGetMatrix(pc, &fact)); 305 PetscCall(MatDuplicate(Bd, MAT_DO_NOT_COPY_VALUES, &AinvBd)); 306 PetscCall(MatMatSolve(fact, Bd, AinvBd)); 307 } else { 308 PetscBool ex = PETSC_TRUE; 309 310 if (ex) { 311 Mat Ainvd; 312 313 PetscCall(PCComputeOperator(pc, MATDENSE, &Ainvd)); 314 PetscCall(MatMatMult(Ainvd, Bd, MAT_INITIAL_MATRIX, fill, &AinvBd)); 315 PetscCall(MatDestroy(&Ainvd)); 316 } else { 317 Vec sol, rhs; 318 PetscScalar *arrayrhs, *arraysol; 319 PetscInt i, nrhs, n; 320 321 PetscCall(MatDuplicate(Bd, MAT_DO_NOT_COPY_VALUES, &AinvBd)); 322 PetscCall(MatGetSize(Bd, &n, &nrhs)); 323 PetscCall(MatDenseGetArray(Bd, &arrayrhs)); 324 PetscCall(MatDenseGetArray(AinvBd, &arraysol)); 325 PetscCall(KSPGetSolution(ksp, &sol)); 326 PetscCall(KSPGetRhs(ksp, &rhs)); 327 for (i = 0; i < nrhs; i++) { 328 PetscCall(VecPlaceArray(rhs, arrayrhs + i * n)); 329 PetscCall(VecPlaceArray(sol, arraysol + i * n)); 330 PetscCall(KSPSolve(ksp, rhs, sol)); 331 PetscCall(VecResetArray(rhs)); 332 PetscCall(VecResetArray(sol)); 333 } 334 PetscCall(MatDenseRestoreArray(Bd, &arrayrhs)); 335 PetscCall(MatDenseRestoreArray(AinvBd, &arrayrhs)); 336 } 337 } 338 if (!Bdense & !issym) PetscCall(MatDestroy(&Bd)); 339 340 if (!issym) { 341 if (!Cdense) { 342 PetscCall(MatConvert(C, MATSEQDENSE, MAT_INITIAL_MATRIX, &Cd)); 343 } else { 344 Cd = C; 345 } 346 PetscCall(MatMatMult(Cd, AinvBd, reuse, fill, S)); 347 if (!Cdense) PetscCall(MatDestroy(&Cd)); 348 } else { 349 PetscCall(MatTransposeMatMult(Bd, AinvBd, reuse, fill, S)); 350 if (!Bdense) PetscCall(MatDestroy(&Bd)); 351 } 352 PetscCall(MatDestroy(&AinvBd)); 353 } 354 355 if (D) { 356 Mat Dd; 357 PetscBool Ddense; 358 359 PetscCall(PetscObjectTypeCompare((PetscObject)D, MATSEQDENSE, &Ddense)); 360 if (!Ddense) { 361 PetscCall(MatConvert(D, MATSEQDENSE, MAT_INITIAL_MATRIX, &Dd)); 362 } else { 363 Dd = D; 364 } 365 if (n_I) { 366 PetscCall(MatAYPX(*S, -1.0, Dd, SAME_NONZERO_PATTERN)); 367 } else { 368 if (reuse == MAT_INITIAL_MATRIX) { 369 PetscCall(MatDuplicate(Dd, MAT_COPY_VALUES, S)); 370 } else { 371 PetscCall(MatCopy(Dd, *S, SAME_NONZERO_PATTERN)); 372 } 373 } 374 if (!Ddense) PetscCall(MatDestroy(&Dd)); 375 } else { 376 PetscCall(MatScale(*S, -1.0)); 377 } 378 PetscFunctionReturn(PETSC_SUCCESS); 379 } 380 381 PetscErrorCode PCBDDCSubSchursSetUp(PCBDDCSubSchurs sub_schurs, Mat Ain, Mat Sin, PetscBool exact_schur, PetscInt xadj[], PetscInt adjncy[], PetscInt nlayers, Vec scaling, PetscBool compute_Stilda, PetscBool reuse_solvers, PetscBool benign_trick, PetscInt benign_n, PetscInt benign_p0_lidx[], IS benign_zerodiag_subs[], Mat change, IS change_primal) 382 { 383 Mat F, A_II, A_IB, A_BI, A_BB, AE_II; 384 Mat S_all; 385 Vec gstash, lstash; 386 VecScatter sstash; 387 IS is_I, is_I_layer; 388 IS all_subsets, all_subsets_mult, all_subsets_n; 389 PetscScalar *stasharray, *Bwork; 390 PetscInt *nnz, *all_local_idx_N; 391 PetscInt *auxnum1, *auxnum2; 392 PetscInt i, subset_size, max_subset_size; 393 PetscInt n_B, extra, local_size, global_size; 394 PetscInt local_stash_size; 395 PetscBLASInt B_N, B_ierr, B_lwork, *pivots; 396 MPI_Comm comm_n; 397 PetscBool deluxe = PETSC_TRUE; 398 PetscBool use_potr = PETSC_FALSE, use_sytr = PETSC_FALSE; 399 PetscViewer matl_dbg_viewer = NULL; 400 PetscBool flg; 401 402 PetscFunctionBegin; 403 PetscCall(MatDestroy(&sub_schurs->A)); 404 PetscCall(MatDestroy(&sub_schurs->S)); 405 if (Ain) { 406 PetscCall(PetscObjectReference((PetscObject)Ain)); 407 sub_schurs->A = Ain; 408 } 409 410 PetscCall(PetscObjectReference((PetscObject)Sin)); 411 sub_schurs->S = Sin; 412 if (sub_schurs->schur_explicit) sub_schurs->schur_explicit = (PetscBool)(!!sub_schurs->A); 413 414 /* preliminary checks */ 415 PetscCheck(sub_schurs->schur_explicit || !compute_Stilda, PetscObjectComm((PetscObject)sub_schurs->l2gmap), PETSC_ERR_SUP, "Adaptive selection of constraints requires MUMPS and/or MKL_PARDISO"); 416 417 if (benign_trick) sub_schurs->is_posdef = PETSC_FALSE; 418 419 /* debug (MATLAB) */ 420 if (sub_schurs->debug) { 421 PetscMPIInt size, rank; 422 PetscInt nr, *print_schurs_ranks, print_schurs = PETSC_FALSE; 423 424 PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)sub_schurs->l2gmap), &size)); 425 PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)sub_schurs->l2gmap), &rank)); 426 nr = size; 427 PetscCall(PetscMalloc1(nr, &print_schurs_ranks)); 428 PetscOptionsBegin(PetscObjectComm((PetscObject)sub_schurs->l2gmap), sub_schurs->prefix, "BDDC sub_schurs options", "PC"); 429 PetscCall(PetscOptionsIntArray("-sub_schurs_debug_ranks", "Ranks to debug (all if the option is not used)", NULL, print_schurs_ranks, &nr, &flg)); 430 if (!flg) print_schurs = PETSC_TRUE; 431 else { 432 print_schurs = PETSC_FALSE; 433 for (i = 0; i < nr; i++) 434 if (print_schurs_ranks[i] == rank) { 435 print_schurs = PETSC_TRUE; 436 break; 437 } 438 } 439 PetscOptionsEnd(); 440 PetscCall(PetscFree(print_schurs_ranks)); 441 if (print_schurs) { 442 char filename[256]; 443 444 PetscCall(PetscSNPrintf(filename, sizeof(filename), "sub_schurs_Schur_r%d.m", PetscGlobalRank)); 445 PetscCall(PetscViewerASCIIOpen(PETSC_COMM_SELF, filename, &matl_dbg_viewer)); 446 PetscCall(PetscViewerPushFormat(matl_dbg_viewer, PETSC_VIEWER_ASCII_MATLAB)); 447 } 448 } 449 450 /* restrict work on active processes */ 451 if (sub_schurs->restrict_comm) { 452 PetscSubcomm subcomm; 453 PetscMPIInt color, rank; 454 455 color = 0; 456 if (!sub_schurs->n_subs) color = 1; /* this can happen if we are in a multilevel case or if the subdomain is disconnected */ 457 PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)sub_schurs->l2gmap), &rank)); 458 PetscCall(PetscSubcommCreate(PetscObjectComm((PetscObject)sub_schurs->l2gmap), &subcomm)); 459 PetscCall(PetscSubcommSetNumber(subcomm, 2)); 460 PetscCall(PetscSubcommSetTypeGeneral(subcomm, color, rank)); 461 PetscCall(PetscCommDuplicate(PetscSubcommChild(subcomm), &comm_n, NULL)); 462 PetscCall(PetscSubcommDestroy(&subcomm)); 463 if (!sub_schurs->n_subs) { 464 PetscCall(PetscCommDestroy(&comm_n)); 465 PetscFunctionReturn(PETSC_SUCCESS); 466 } 467 } else { 468 PetscCall(PetscCommDuplicate(PetscObjectComm((PetscObject)sub_schurs->l2gmap), &comm_n, NULL)); 469 } 470 471 /* get Schur complement matrices */ 472 if (!sub_schurs->schur_explicit) { 473 Mat tA_IB, tA_BI, tA_BB; 474 PetscBool isseqsbaij; 475 PetscCall(MatSchurComplementGetSubMatrices(sub_schurs->S, &A_II, NULL, &tA_IB, &tA_BI, &tA_BB)); 476 PetscCall(PetscObjectTypeCompare((PetscObject)tA_BB, MATSEQSBAIJ, &isseqsbaij)); 477 if (isseqsbaij) { 478 PetscCall(MatConvert(tA_BB, MATSEQAIJ, MAT_INITIAL_MATRIX, &A_BB)); 479 PetscCall(MatConvert(tA_IB, MATSEQAIJ, MAT_INITIAL_MATRIX, &A_IB)); 480 PetscCall(MatConvert(tA_BI, MATSEQAIJ, MAT_INITIAL_MATRIX, &A_BI)); 481 } else { 482 PetscCall(PetscObjectReference((PetscObject)tA_BB)); 483 A_BB = tA_BB; 484 PetscCall(PetscObjectReference((PetscObject)tA_IB)); 485 A_IB = tA_IB; 486 PetscCall(PetscObjectReference((PetscObject)tA_BI)); 487 A_BI = tA_BI; 488 } 489 } else { 490 A_II = NULL; 491 A_IB = NULL; 492 A_BI = NULL; 493 A_BB = NULL; 494 } 495 S_all = NULL; 496 497 /* determine interior problems */ 498 PetscCall(ISGetLocalSize(sub_schurs->is_I, &i)); 499 if (nlayers >= 0 && i) { /* Interior problems can be different from the original one */ 500 PetscBT touched; 501 const PetscInt *idx_B; 502 PetscInt n_I, n_B, n_local_dofs, n_prev_added, j, layer, *local_numbering; 503 504 PetscCheck(xadj, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Cannot request layering without adjacency"); 505 /* get sizes */ 506 PetscCall(ISGetLocalSize(sub_schurs->is_I, &n_I)); 507 PetscCall(ISGetLocalSize(sub_schurs->is_B, &n_B)); 508 509 PetscCall(PetscMalloc1(n_I + n_B, &local_numbering)); 510 PetscCall(PetscBTCreate(n_I + n_B, &touched)); 511 PetscCall(PetscBTMemzero(n_I + n_B, touched)); 512 513 /* all boundary dofs must be skipped when adding layers */ 514 PetscCall(ISGetIndices(sub_schurs->is_B, &idx_B)); 515 for (j = 0; j < n_B; j++) PetscCall(PetscBTSet(touched, idx_B[j])); 516 PetscCall(PetscArraycpy(local_numbering, idx_B, n_B)); 517 PetscCall(ISRestoreIndices(sub_schurs->is_B, &idx_B)); 518 519 /* add prescribed number of layers of dofs */ 520 n_local_dofs = n_B; 521 n_prev_added = n_B; 522 for (layer = 0; layer < nlayers; layer++) { 523 PetscInt n_added = 0; 524 if (n_local_dofs == n_I + n_B) break; 525 PetscCheck(n_local_dofs <= n_I + n_B, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Error querying layer %" PetscInt_FMT ". Out of bound access (%" PetscInt_FMT " > %" PetscInt_FMT ")", layer, n_local_dofs, n_I + n_B); 526 PetscCall(PCBDDCAdjGetNextLayer_Private(local_numbering + n_local_dofs, n_prev_added, touched, xadj, adjncy, &n_added)); 527 n_prev_added = n_added; 528 n_local_dofs += n_added; 529 if (!n_added) break; 530 } 531 PetscCall(PetscBTDestroy(&touched)); 532 533 /* IS for I layer dofs in original numbering */ 534 PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)sub_schurs->is_I), n_local_dofs - n_B, local_numbering + n_B, PETSC_COPY_VALUES, &is_I_layer)); 535 PetscCall(PetscFree(local_numbering)); 536 PetscCall(ISSort(is_I_layer)); 537 /* IS for I layer dofs in I numbering */ 538 if (!sub_schurs->schur_explicit) { 539 ISLocalToGlobalMapping ItoNmap; 540 PetscCall(ISLocalToGlobalMappingCreateIS(sub_schurs->is_I, &ItoNmap)); 541 PetscCall(ISGlobalToLocalMappingApplyIS(ItoNmap, IS_GTOLM_DROP, is_I_layer, &is_I)); 542 PetscCall(ISLocalToGlobalMappingDestroy(&ItoNmap)); 543 544 /* II block */ 545 PetscCall(MatCreateSubMatrix(A_II, is_I, is_I, MAT_INITIAL_MATRIX, &AE_II)); 546 } 547 } else { 548 PetscInt n_I; 549 550 /* IS for I dofs in original numbering */ 551 PetscCall(PetscObjectReference((PetscObject)sub_schurs->is_I)); 552 is_I_layer = sub_schurs->is_I; 553 554 /* IS for I dofs in I numbering (strided 1) */ 555 if (!sub_schurs->schur_explicit) { 556 PetscCall(ISGetSize(sub_schurs->is_I, &n_I)); 557 PetscCall(ISCreateStride(PetscObjectComm((PetscObject)sub_schurs->is_I), n_I, 0, 1, &is_I)); 558 559 /* II block is the same */ 560 PetscCall(PetscObjectReference((PetscObject)A_II)); 561 AE_II = A_II; 562 } 563 } 564 565 /* Get info on subset sizes and sum of all subsets sizes */ 566 max_subset_size = 0; 567 local_size = 0; 568 for (i = 0; i < sub_schurs->n_subs; i++) { 569 PetscCall(ISGetLocalSize(sub_schurs->is_subs[i], &subset_size)); 570 max_subset_size = PetscMax(subset_size, max_subset_size); 571 local_size += subset_size; 572 } 573 574 /* Work arrays for local indices */ 575 extra = 0; 576 PetscCall(ISGetLocalSize(sub_schurs->is_B, &n_B)); 577 if (sub_schurs->schur_explicit && is_I_layer) PetscCall(ISGetLocalSize(is_I_layer, &extra)); 578 PetscCall(PetscMalloc1(n_B + extra, &all_local_idx_N)); 579 if (extra) { 580 const PetscInt *idxs; 581 PetscCall(ISGetIndices(is_I_layer, &idxs)); 582 PetscCall(PetscArraycpy(all_local_idx_N, idxs, extra)); 583 PetscCall(ISRestoreIndices(is_I_layer, &idxs)); 584 } 585 PetscCall(PetscMalloc1(sub_schurs->n_subs, &auxnum1)); 586 PetscCall(PetscMalloc1(sub_schurs->n_subs, &auxnum2)); 587 588 /* Get local indices in local numbering */ 589 local_size = 0; 590 local_stash_size = 0; 591 for (i = 0; i < sub_schurs->n_subs; i++) { 592 const PetscInt *idxs; 593 594 PetscCall(ISGetLocalSize(sub_schurs->is_subs[i], &subset_size)); 595 PetscCall(ISGetIndices(sub_schurs->is_subs[i], &idxs)); 596 /* start (smallest in global ordering) and multiplicity */ 597 auxnum1[i] = idxs[0]; 598 auxnum2[i] = subset_size * subset_size; 599 /* subset indices in local numbering */ 600 PetscCall(PetscArraycpy(all_local_idx_N + local_size + extra, idxs, subset_size)); 601 PetscCall(ISRestoreIndices(sub_schurs->is_subs[i], &idxs)); 602 local_size += subset_size; 603 local_stash_size += subset_size * subset_size; 604 } 605 606 /* allocate extra workspace needed only for GETRI or SYTRF when inverting the blocks or the entire Schur complement */ 607 use_potr = use_sytr = PETSC_FALSE; 608 if (benign_trick || (sub_schurs->is_hermitian && sub_schurs->is_posdef)) { 609 use_potr = PETSC_TRUE; 610 } else if (sub_schurs->is_symmetric) { 611 use_sytr = PETSC_TRUE; 612 } 613 if (local_size && !use_potr && compute_Stilda) { 614 PetscScalar lwork, dummyscalar = 0.; 615 PetscBLASInt dummyint = 0; 616 617 B_lwork = -1; 618 PetscCall(PetscBLASIntCast(local_size, &B_N)); 619 PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); 620 if (use_sytr) { 621 PetscCallBLAS("LAPACKsytrf", LAPACKsytrf_("L", &B_N, &dummyscalar, &B_N, &dummyint, &lwork, &B_lwork, &B_ierr)); 622 PetscCheck(!B_ierr, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in query to SYTRF Lapack routine %" PetscBLASInt_FMT, B_ierr); 623 } else { 624 PetscCallBLAS("LAPACKgetri", LAPACKgetri_(&B_N, &dummyscalar, &B_N, &dummyint, &lwork, &B_lwork, &B_ierr)); 625 PetscCheck(!B_ierr, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in query to GETRI Lapack routine %" PetscBLASInt_FMT, B_ierr); 626 } 627 PetscCall(PetscFPTrapPop()); 628 PetscCall(PetscBLASIntCast((PetscInt)PetscRealPart(lwork), &B_lwork)); 629 PetscCall(PetscMalloc2(B_lwork, &Bwork, B_N, &pivots)); 630 } else { 631 Bwork = NULL; 632 pivots = NULL; 633 } 634 635 /* prepare data for summing up properly schurs on subsets */ 636 PetscCall(ISCreateGeneral(comm_n, sub_schurs->n_subs, auxnum1, PETSC_OWN_POINTER, &all_subsets_n)); 637 PetscCall(ISLocalToGlobalMappingApplyIS(sub_schurs->l2gmap, all_subsets_n, &all_subsets)); 638 PetscCall(ISDestroy(&all_subsets_n)); 639 PetscCall(ISCreateGeneral(comm_n, sub_schurs->n_subs, auxnum2, PETSC_OWN_POINTER, &all_subsets_mult)); 640 PetscCall(ISRenumber(all_subsets, all_subsets_mult, &global_size, &all_subsets_n)); 641 PetscCall(ISDestroy(&all_subsets)); 642 PetscCall(ISDestroy(&all_subsets_mult)); 643 PetscCall(ISGetLocalSize(all_subsets_n, &i)); 644 PetscCheck(i == local_stash_size, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid size of new subset! %" PetscInt_FMT " != %" PetscInt_FMT, i, local_stash_size); 645 PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, local_stash_size, NULL, &lstash)); 646 PetscCall(VecCreateMPI(comm_n, PETSC_DECIDE, global_size, &gstash)); 647 PetscCall(VecScatterCreate(lstash, NULL, gstash, all_subsets_n, &sstash)); 648 PetscCall(ISDestroy(&all_subsets_n)); 649 650 /* subset indices in local boundary numbering */ 651 if (!sub_schurs->is_Ej_all) { 652 PetscInt *all_local_idx_B; 653 654 PetscCall(PetscMalloc1(local_size, &all_local_idx_B)); 655 PetscCall(ISGlobalToLocalMappingApply(sub_schurs->BtoNmap, IS_GTOLM_DROP, local_size, all_local_idx_N + extra, &subset_size, all_local_idx_B)); 656 PetscCheck(subset_size == local_size, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Error in sub_schurs serial (BtoNmap)! %" PetscInt_FMT " != %" PetscInt_FMT, subset_size, local_size); 657 PetscCall(ISCreateGeneral(PETSC_COMM_SELF, local_size, all_local_idx_B, PETSC_OWN_POINTER, &sub_schurs->is_Ej_all)); 658 } 659 660 if (change) { 661 ISLocalToGlobalMapping BtoS; 662 IS change_primal_B; 663 IS change_primal_all; 664 665 PetscCheck(!sub_schurs->change_primal_sub, PETSC_COMM_SELF, PETSC_ERR_PLIB, "This should not happen"); 666 PetscCheck(!sub_schurs->change, PETSC_COMM_SELF, PETSC_ERR_PLIB, "This should not happen"); 667 PetscCall(PetscMalloc1(sub_schurs->n_subs, &sub_schurs->change_primal_sub)); 668 for (i = 0; i < sub_schurs->n_subs; i++) { 669 ISLocalToGlobalMapping NtoS; 670 PetscCall(ISLocalToGlobalMappingCreateIS(sub_schurs->is_subs[i], &NtoS)); 671 PetscCall(ISGlobalToLocalMappingApplyIS(NtoS, IS_GTOLM_DROP, change_primal, &sub_schurs->change_primal_sub[i])); 672 PetscCall(ISLocalToGlobalMappingDestroy(&NtoS)); 673 } 674 PetscCall(ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap, IS_GTOLM_DROP, change_primal, &change_primal_B)); 675 PetscCall(ISLocalToGlobalMappingCreateIS(sub_schurs->is_Ej_all, &BtoS)); 676 PetscCall(ISGlobalToLocalMappingApplyIS(BtoS, IS_GTOLM_DROP, change_primal_B, &change_primal_all)); 677 PetscCall(ISLocalToGlobalMappingDestroy(&BtoS)); 678 PetscCall(ISDestroy(&change_primal_B)); 679 PetscCall(PetscMalloc1(sub_schurs->n_subs, &sub_schurs->change)); 680 for (i = 0; i < sub_schurs->n_subs; i++) { 681 Mat change_sub; 682 683 PetscCall(ISGetLocalSize(sub_schurs->is_subs[i], &subset_size)); 684 PetscCall(KSPCreate(PETSC_COMM_SELF, &sub_schurs->change[i])); 685 PetscCall(KSPSetNestLevel(sub_schurs->change[i], 1)); /* do not seem to have direct access to a PC from which to get the level of nests */ 686 PetscCall(KSPSetType(sub_schurs->change[i], KSPPREONLY)); 687 if (!sub_schurs->change_with_qr) { 688 PetscCall(MatCreateSubMatrix(change, sub_schurs->is_subs[i], sub_schurs->is_subs[i], MAT_INITIAL_MATRIX, &change_sub)); 689 } else { 690 Mat change_subt; 691 PetscCall(MatCreateSubMatrix(change, sub_schurs->is_subs[i], sub_schurs->is_subs[i], MAT_INITIAL_MATRIX, &change_subt)); 692 PetscCall(MatConvert(change_subt, MATSEQDENSE, MAT_INITIAL_MATRIX, &change_sub)); 693 PetscCall(MatDestroy(&change_subt)); 694 } 695 PetscCall(KSPSetOperators(sub_schurs->change[i], change_sub, change_sub)); 696 PetscCall(MatDestroy(&change_sub)); 697 PetscCall(KSPSetOptionsPrefix(sub_schurs->change[i], sub_schurs->prefix)); 698 PetscCall(KSPAppendOptionsPrefix(sub_schurs->change[i], "sub_schurs_change_")); 699 } 700 PetscCall(ISDestroy(&change_primal_all)); 701 } 702 703 /* Local matrix of all local Schur on subsets (transposed) */ 704 if (!sub_schurs->S_Ej_all) { 705 Mat T; 706 PetscScalar *v; 707 PetscInt *ii, *jj; 708 PetscInt cum, i, j, k; 709 710 /* MatSeqAIJSetPreallocation + MatSetValues is slow for these kind of matrices (may have large blocks) 711 Allocate properly a representative matrix and duplicate */ 712 PetscCall(PetscMalloc3(local_size + 1, &ii, local_stash_size, &jj, local_stash_size, &v)); 713 ii[0] = 0; 714 cum = 0; 715 for (i = 0; i < sub_schurs->n_subs; i++) { 716 PetscCall(ISGetLocalSize(sub_schurs->is_subs[i], &subset_size)); 717 for (j = 0; j < subset_size; j++) { 718 const PetscInt row = cum + j; 719 PetscInt col = cum; 720 721 ii[row + 1] = ii[row] + subset_size; 722 for (k = ii[row]; k < ii[row + 1]; k++) { 723 jj[k] = col; 724 col++; 725 } 726 } 727 cum += subset_size; 728 } 729 PetscCall(MatCreateSeqAIJWithArrays(PETSC_COMM_SELF, local_size, local_size, ii, jj, v, &T)); 730 PetscCall(MatDuplicate(T, MAT_DO_NOT_COPY_VALUES, &sub_schurs->S_Ej_all)); 731 PetscCall(MatDestroy(&T)); 732 PetscCall(PetscFree3(ii, jj, v)); 733 } 734 /* matrices for deluxe scaling and adaptive selection */ 735 if (compute_Stilda) { 736 if (!sub_schurs->sum_S_Ej_tilda_all) PetscCall(MatDuplicate(sub_schurs->S_Ej_all, MAT_DO_NOT_COPY_VALUES, &sub_schurs->sum_S_Ej_tilda_all)); 737 if (!sub_schurs->sum_S_Ej_inv_all && deluxe) PetscCall(MatDuplicate(sub_schurs->S_Ej_all, MAT_DO_NOT_COPY_VALUES, &sub_schurs->sum_S_Ej_inv_all)); 738 } 739 740 /* Compute Schur complements explicitly */ 741 F = NULL; 742 if (!sub_schurs->schur_explicit) { 743 /* this code branch is used when MatFactor with Schur complement support is not present or when explicitly requested; 744 it is not efficient, unless the economic version of the scaling is used */ 745 Mat S_Ej_expl; 746 PetscScalar *work; 747 PetscInt j, *dummy_idx; 748 PetscBool Sdense; 749 750 PetscCall(PetscMalloc2(max_subset_size, &dummy_idx, max_subset_size * max_subset_size, &work)); 751 local_size = 0; 752 for (i = 0; i < sub_schurs->n_subs; i++) { 753 IS is_subset_B; 754 Mat AE_EE, AE_IE, AE_EI, S_Ej; 755 756 /* subsets in original and boundary numbering */ 757 PetscCall(ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap, IS_GTOLM_DROP, sub_schurs->is_subs[i], &is_subset_B)); 758 /* EE block */ 759 PetscCall(MatCreateSubMatrix(A_BB, is_subset_B, is_subset_B, MAT_INITIAL_MATRIX, &AE_EE)); 760 /* IE block */ 761 PetscCall(MatCreateSubMatrix(A_IB, is_I, is_subset_B, MAT_INITIAL_MATRIX, &AE_IE)); 762 /* EI block */ 763 if (sub_schurs->is_symmetric) { 764 PetscCall(MatCreateTranspose(AE_IE, &AE_EI)); 765 } else if (sub_schurs->is_hermitian) { 766 PetscCall(MatCreateHermitianTranspose(AE_IE, &AE_EI)); 767 } else { 768 PetscCall(MatCreateSubMatrix(A_BI, is_subset_B, is_I, MAT_INITIAL_MATRIX, &AE_EI)); 769 } 770 PetscCall(ISDestroy(&is_subset_B)); 771 PetscCall(MatCreateSchurComplement(AE_II, AE_II, AE_IE, AE_EI, AE_EE, &S_Ej)); 772 PetscCall(MatDestroy(&AE_EE)); 773 PetscCall(MatDestroy(&AE_IE)); 774 PetscCall(MatDestroy(&AE_EI)); 775 if (AE_II == A_II) { /* we can reuse the same ksp */ 776 KSP ksp; 777 PetscCall(MatSchurComplementGetKSP(sub_schurs->S, &ksp)); 778 PetscCall(MatSchurComplementSetKSP(S_Ej, ksp)); 779 } else { /* build new ksp object which inherits ksp and pc types from the original one */ 780 KSP origksp, schurksp; 781 PC origpc, schurpc; 782 KSPType ksp_type; 783 PetscInt n_internal; 784 PetscBool ispcnone; 785 786 PetscCall(MatSchurComplementGetKSP(sub_schurs->S, &origksp)); 787 PetscCall(MatSchurComplementGetKSP(S_Ej, &schurksp)); 788 PetscCall(KSPGetType(origksp, &ksp_type)); 789 PetscCall(KSPSetType(schurksp, ksp_type)); 790 PetscCall(KSPGetPC(schurksp, &schurpc)); 791 PetscCall(KSPGetPC(origksp, &origpc)); 792 PetscCall(PetscObjectTypeCompare((PetscObject)origpc, PCNONE, &ispcnone)); 793 if (!ispcnone) { 794 PCType pc_type; 795 PetscCall(PCGetType(origpc, &pc_type)); 796 PetscCall(PCSetType(schurpc, pc_type)); 797 } else { 798 PetscCall(PCSetType(schurpc, PCLU)); 799 } 800 PetscCall(ISGetSize(is_I, &n_internal)); 801 if (!n_internal) { /* UMFPACK gives error with 0 sized problems */ 802 MatSolverType solver = NULL; 803 PetscCall(PCFactorGetMatSolverType(origpc, &solver)); 804 if (solver) PetscCall(PCFactorSetMatSolverType(schurpc, solver)); 805 } 806 PetscCall(KSPSetUp(schurksp)); 807 } 808 PetscCall(ISGetLocalSize(sub_schurs->is_subs[i], &subset_size)); 809 PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, subset_size, subset_size, work, &S_Ej_expl)); 810 PetscCall(PCBDDCComputeExplicitSchur(S_Ej, sub_schurs->is_symmetric, MAT_REUSE_MATRIX, &S_Ej_expl)); 811 PetscCall(PetscObjectTypeCompare((PetscObject)S_Ej_expl, MATSEQDENSE, &Sdense)); 812 PetscCheck(Sdense, PETSC_COMM_SELF, PETSC_ERR_SUP, "Not yet implemented for sparse matrices"); 813 for (j = 0; j < subset_size; j++) dummy_idx[j] = local_size + j; 814 PetscCall(MatSetValues(sub_schurs->S_Ej_all, subset_size, dummy_idx, subset_size, dummy_idx, work, INSERT_VALUES)); 815 PetscCall(MatDestroy(&S_Ej)); 816 PetscCall(MatDestroy(&S_Ej_expl)); 817 local_size += subset_size; 818 } 819 PetscCall(PetscFree2(dummy_idx, work)); 820 /* free */ 821 PetscCall(ISDestroy(&is_I)); 822 PetscCall(MatDestroy(&AE_II)); 823 PetscCall(PetscFree(all_local_idx_N)); 824 } else { 825 Mat A, cs_AIB_mat = NULL, benign_AIIm1_ones_mat = NULL; 826 Mat *gdswA; 827 Vec Dall = NULL; 828 IS is_A_all, *is_p_r = NULL; 829 MatType Stype; 830 PetscScalar *work, *S_data, *schur_factor, infty = PETSC_MAX_REAL; 831 PetscScalar *SEj_arr = NULL, *SEjinv_arr = NULL; 832 const PetscScalar *rS_data; 833 PetscInt n, n_I, size_schur, size_active_schur, cum, cum2; 834 PetscBool economic, solver_S, S_lower_triangular = PETSC_FALSE; 835 PetscBool schur_has_vertices, factor_workaround; 836 PetscBool use_cholesky; 837 #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA) 838 PetscBool oldpin; 839 #endif 840 841 /* get sizes */ 842 n_I = 0; 843 if (is_I_layer) PetscCall(ISGetLocalSize(is_I_layer, &n_I)); 844 economic = PETSC_FALSE; 845 PetscCall(ISGetLocalSize(sub_schurs->is_I, &cum)); 846 if (cum != n_I) economic = PETSC_TRUE; 847 PetscCall(MatGetLocalSize(sub_schurs->A, &n, NULL)); 848 size_active_schur = local_size; 849 850 /* import scaling vector (wrong formulation if we have 3D edges) */ 851 if (scaling && compute_Stilda) { 852 const PetscScalar *array; 853 PetscScalar *array2; 854 const PetscInt *idxs; 855 PetscInt i; 856 857 PetscCall(ISGetIndices(sub_schurs->is_Ej_all, &idxs)); 858 PetscCall(VecCreateSeq(PETSC_COMM_SELF, size_active_schur, &Dall)); 859 PetscCall(VecGetArrayRead(scaling, &array)); 860 PetscCall(VecGetArray(Dall, &array2)); 861 for (i = 0; i < size_active_schur; i++) array2[i] = array[idxs[i]]; 862 PetscCall(VecRestoreArray(Dall, &array2)); 863 PetscCall(VecRestoreArrayRead(scaling, &array)); 864 PetscCall(ISRestoreIndices(sub_schurs->is_Ej_all, &idxs)); 865 deluxe = PETSC_FALSE; 866 } 867 868 /* size active schurs does not count any dirichlet or vertex dof on the interface */ 869 factor_workaround = PETSC_FALSE; 870 schur_has_vertices = PETSC_FALSE; 871 cum = n_I + size_active_schur; 872 if (sub_schurs->is_dir) { 873 const PetscInt *idxs; 874 PetscInt n_dir; 875 876 PetscCall(ISGetLocalSize(sub_schurs->is_dir, &n_dir)); 877 PetscCall(ISGetIndices(sub_schurs->is_dir, &idxs)); 878 PetscCall(PetscArraycpy(all_local_idx_N + cum, idxs, n_dir)); 879 PetscCall(ISRestoreIndices(sub_schurs->is_dir, &idxs)); 880 cum += n_dir; 881 if (!sub_schurs->gdsw) factor_workaround = PETSC_TRUE; 882 } 883 /* include the primal vertices in the Schur complement */ 884 if (exact_schur && sub_schurs->is_vertices && (compute_Stilda || benign_n)) { 885 PetscInt n_v; 886 887 PetscCall(ISGetLocalSize(sub_schurs->is_vertices, &n_v)); 888 if (n_v) { 889 const PetscInt *idxs; 890 891 PetscCall(ISGetIndices(sub_schurs->is_vertices, &idxs)); 892 PetscCall(PetscArraycpy(all_local_idx_N + cum, idxs, n_v)); 893 PetscCall(ISRestoreIndices(sub_schurs->is_vertices, &idxs)); 894 cum += n_v; 895 if (!sub_schurs->gdsw) factor_workaround = PETSC_TRUE; 896 schur_has_vertices = PETSC_TRUE; 897 } 898 } 899 size_schur = cum - n_I; 900 PetscCall(ISCreateGeneral(PETSC_COMM_SELF, cum, all_local_idx_N, PETSC_OWN_POINTER, &is_A_all)); 901 #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA) 902 oldpin = sub_schurs->A->boundtocpu; 903 PetscCall(MatBindToCPU(sub_schurs->A, PETSC_TRUE)); 904 #endif 905 if (cum == n) { 906 PetscCall(ISSetPermutation(is_A_all)); 907 PetscCall(MatPermute(sub_schurs->A, is_A_all, is_A_all, &A)); 908 } else { 909 PetscCall(MatCreateSubMatrix(sub_schurs->A, is_A_all, is_A_all, MAT_INITIAL_MATRIX, &A)); 910 } 911 #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA) 912 PetscCall(MatBindToCPU(sub_schurs->A, oldpin)); 913 #endif 914 PetscCall(MatSetOptionsPrefixFactor(A, sub_schurs->prefix)); 915 PetscCall(MatAppendOptionsPrefixFactor(A, "sub_schurs_")); 916 917 /* if we actually change the basis for the pressures, LDL^T factors will use a lot of memory 918 this is a workaround */ 919 if (benign_n) { 920 Vec v, benign_AIIm1_ones; 921 ISLocalToGlobalMapping N_to_reor; 922 IS is_p0, is_p0_p; 923 PetscScalar *cs_AIB, *AIIm1_data; 924 PetscInt sizeA; 925 926 PetscCall(ISLocalToGlobalMappingCreateIS(is_A_all, &N_to_reor)); 927 PetscCall(ISCreateGeneral(PETSC_COMM_SELF, benign_n, benign_p0_lidx, PETSC_COPY_VALUES, &is_p0)); 928 PetscCall(ISGlobalToLocalMappingApplyIS(N_to_reor, IS_GTOLM_DROP, is_p0, &is_p0_p)); 929 PetscCall(ISDestroy(&is_p0)); 930 PetscCall(MatCreateVecs(A, &v, &benign_AIIm1_ones)); 931 PetscCall(VecGetSize(v, &sizeA)); 932 PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, sizeA, benign_n, NULL, &benign_AIIm1_ones_mat)); 933 PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, size_schur, benign_n, NULL, &cs_AIB_mat)); 934 PetscCall(MatDenseGetArray(cs_AIB_mat, &cs_AIB)); 935 PetscCall(MatDenseGetArray(benign_AIIm1_ones_mat, &AIIm1_data)); 936 PetscCall(PetscMalloc1(benign_n, &is_p_r)); 937 /* compute colsum of A_IB restricted to pressures */ 938 for (i = 0; i < benign_n; i++) { 939 const PetscScalar *array; 940 const PetscInt *idxs; 941 PetscInt j, nz; 942 943 PetscCall(ISGlobalToLocalMappingApplyIS(N_to_reor, IS_GTOLM_DROP, benign_zerodiag_subs[i], &is_p_r[i])); 944 PetscCall(ISGetLocalSize(is_p_r[i], &nz)); 945 PetscCall(ISGetIndices(is_p_r[i], &idxs)); 946 for (j = 0; j < nz; j++) AIIm1_data[idxs[j] + sizeA * i] = 1.; 947 PetscCall(ISRestoreIndices(is_p_r[i], &idxs)); 948 PetscCall(VecPlaceArray(benign_AIIm1_ones, AIIm1_data + sizeA * i)); 949 PetscCall(MatMult(A, benign_AIIm1_ones, v)); 950 PetscCall(VecResetArray(benign_AIIm1_ones)); 951 PetscCall(VecGetArrayRead(v, &array)); 952 for (j = 0; j < size_schur; j++) { 953 #if defined(PETSC_USE_COMPLEX) 954 cs_AIB[i * size_schur + j] = (PetscRealPart(array[j + n_I]) / nz + PETSC_i * (PetscImaginaryPart(array[j + n_I]) / nz)); 955 #else 956 cs_AIB[i * size_schur + j] = array[j + n_I] / nz; 957 #endif 958 } 959 PetscCall(VecRestoreArrayRead(v, &array)); 960 } 961 PetscCall(MatDenseRestoreArray(cs_AIB_mat, &cs_AIB)); 962 PetscCall(MatDenseRestoreArray(benign_AIIm1_ones_mat, &AIIm1_data)); 963 PetscCall(VecDestroy(&v)); 964 PetscCall(VecDestroy(&benign_AIIm1_ones)); 965 PetscCall(MatSetOption(A, MAT_KEEP_NONZERO_PATTERN, PETSC_FALSE)); 966 PetscCall(MatSetOption(A, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_FALSE)); 967 PetscCall(MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE)); 968 PetscCall(MatZeroRowsColumnsIS(A, is_p0_p, 1.0, NULL, NULL)); 969 PetscCall(ISDestroy(&is_p0_p)); 970 PetscCall(ISLocalToGlobalMappingDestroy(&N_to_reor)); 971 } 972 PetscCall(MatSetOption(A, MAT_SYMMETRIC, sub_schurs->is_symmetric)); 973 PetscCall(MatSetOption(A, MAT_HERMITIAN, sub_schurs->is_hermitian)); 974 PetscCall(MatSetOption(A, MAT_SPD, sub_schurs->is_posdef)); 975 976 /* for complexes, symmetric and hermitian at the same time implies null imaginary part */ 977 use_cholesky = (PetscBool)((use_potr || use_sytr) && sub_schurs->is_hermitian && sub_schurs->is_symmetric); 978 979 /* when using the benign subspace trick, the local Schur complements are SPD */ 980 /* MKL_PARDISO does not handle well the computation of a Schur complement from a symmetric indefinite factorization 981 Use LU and adapt pivoting perturbation (still, solution is not as accurate as with using MUMPS) */ 982 if (benign_trick) { 983 sub_schurs->is_posdef = PETSC_TRUE; 984 PetscCall(PetscStrcmp(sub_schurs->mat_solver_type, MATSOLVERMKL_PARDISO, &flg)); 985 if (flg) use_cholesky = PETSC_FALSE; 986 } 987 if (sub_schurs->mat_factor_type == MAT_FACTOR_NONE) sub_schurs->mat_factor_type = use_cholesky ? MAT_FACTOR_CHOLESKY : MAT_FACTOR_LU; 988 989 if (n_I) { 990 IS is_schur; 991 char stype[64]; 992 PetscBool gpu = PETSC_FALSE; 993 994 PetscCall(MatGetFactor(A, sub_schurs->mat_solver_type, sub_schurs->mat_factor_type, &F)); 995 PetscCheck(F, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "MatGetFactor not supported by matrix instance of type %s. Rerun with \"-info :mat | grep MatGetFactor_\" for additional information", ((PetscObject)A)->type_name); 996 PetscCall(MatSetErrorIfFailure(A, PETSC_TRUE)); 997 #if defined(PETSC_HAVE_MKL_PARDISO) 998 if (benign_trick) PetscCall(MatMkl_PardisoSetCntl(F, 10, 10)); 999 #endif 1000 /* subsets ordered last */ 1001 PetscCall(ISCreateStride(PETSC_COMM_SELF, size_schur, n_I, 1, &is_schur)); 1002 PetscCall(MatFactorSetSchurIS(F, is_schur)); 1003 PetscCall(ISDestroy(&is_schur)); 1004 1005 /* factorization step */ 1006 switch (sub_schurs->mat_factor_type) { 1007 case MAT_FACTOR_CHOLESKY: 1008 PetscCall(MatCholeskyFactorSymbolic(F, A, NULL, NULL)); 1009 /* be sure that icntl 19 is not set by command line */ 1010 PetscCall(MatMumpsSetIcntl(F, 19, 2)); 1011 PetscCall(MatCholeskyFactorNumeric(F, A, NULL)); 1012 S_lower_triangular = PETSC_TRUE; 1013 break; 1014 case MAT_FACTOR_LU: 1015 PetscCall(MatLUFactorSymbolic(F, A, NULL, NULL, NULL)); 1016 /* be sure that icntl 19 is not set by command line */ 1017 PetscCall(MatMumpsSetIcntl(F, 19, 3)); 1018 PetscCall(MatLUFactorNumeric(F, A, NULL)); 1019 break; 1020 default: 1021 SETERRQ(PetscObjectComm((PetscObject)F), PETSC_ERR_SUP, "Unsupported factor type %s", MatFactorTypes[sub_schurs->mat_factor_type]); 1022 } 1023 PetscCall(MatViewFromOptions(F, (PetscObject)A, "-mat_factor_view")); 1024 1025 if (matl_dbg_viewer) { 1026 Mat S; 1027 IS is; 1028 1029 PetscCall(PetscObjectSetName((PetscObject)A, "A")); 1030 PetscCall(MatView(A, matl_dbg_viewer)); 1031 PetscCall(MatFactorCreateSchurComplement(F, &S, NULL)); 1032 PetscCall(PetscObjectSetName((PetscObject)S, "S")); 1033 PetscCall(MatView(S, matl_dbg_viewer)); 1034 PetscCall(MatDestroy(&S)); 1035 PetscCall(ISCreateStride(PETSC_COMM_SELF, n_I, 0, 1, &is)); 1036 PetscCall(PetscObjectSetName((PetscObject)is, "I")); 1037 PetscCall(ISView(is, matl_dbg_viewer)); 1038 PetscCall(ISDestroy(&is)); 1039 PetscCall(ISCreateStride(PETSC_COMM_SELF, size_schur, n_I, 1, &is)); 1040 PetscCall(PetscObjectSetName((PetscObject)is, "B")); 1041 PetscCall(ISView(is, matl_dbg_viewer)); 1042 PetscCall(ISDestroy(&is)); 1043 PetscCall(PetscObjectSetName((PetscObject)is_A_all, "IA")); 1044 PetscCall(ISView(is_A_all, matl_dbg_viewer)); 1045 for (i = 0, cum = 0; i < sub_schurs->n_subs; i++) { 1046 IS is; 1047 char name[16]; 1048 1049 PetscCall(PetscSNPrintf(name, sizeof(name), "IE%" PetscInt_FMT, i)); 1050 PetscCall(ISGetLocalSize(sub_schurs->is_subs[i], &subset_size)); 1051 PetscCall(ISCreateStride(PETSC_COMM_SELF, subset_size, cum, 1, &is)); 1052 PetscCall(PetscObjectSetName((PetscObject)is, name)); 1053 PetscCall(ISView(is, matl_dbg_viewer)); 1054 PetscCall(ISDestroy(&is)); 1055 if (sub_schurs->change) { 1056 Mat T; 1057 1058 PetscCall(PetscSNPrintf(name, sizeof(name), "TE%" PetscInt_FMT, i)); 1059 PetscCall(KSPGetOperators(sub_schurs->change[i], &T, NULL)); 1060 PetscCall(PetscObjectSetName((PetscObject)T, name)); 1061 PetscCall(MatView(T, matl_dbg_viewer)); 1062 PetscCall(PetscSNPrintf(name, sizeof(name), "ITE%" PetscInt_FMT, i)); 1063 PetscCall(PetscObjectSetName((PetscObject)sub_schurs->change_primal_sub[i], name)); 1064 PetscCall(ISView(sub_schurs->change_primal_sub[i], matl_dbg_viewer)); 1065 } 1066 cum += subset_size; 1067 } 1068 PetscCall(PetscViewerFlush(matl_dbg_viewer)); 1069 } 1070 1071 /* get explicit Schur Complement computed during numeric factorization */ 1072 PetscCall(MatFactorGetSchurComplement(F, &S_all, NULL)); 1073 PetscCall(PetscStrncpy(stype, MATSEQDENSE, sizeof(stype))); 1074 #if defined(PETSC_HAVE_CUDA) 1075 PetscCall(PetscObjectTypeCompareAny((PetscObject)A, &gpu, MATSEQAIJVIENNACL, MATSEQAIJCUSPARSE, "")); 1076 #endif 1077 if (gpu) PetscCall(PetscStrncpy(stype, MATSEQDENSECUDA, sizeof(stype))); 1078 PetscCall(PetscOptionsGetString(NULL, sub_schurs->prefix, "-sub_schurs_schur_mat_type", stype, sizeof(stype), NULL)); 1079 PetscCall(MatConvert(S_all, stype, MAT_INPLACE_MATRIX, &S_all)); 1080 PetscCall(MatSetOption(S_all, MAT_SPD, sub_schurs->is_posdef)); 1081 PetscCall(MatSetOption(S_all, MAT_HERMITIAN, sub_schurs->is_hermitian)); 1082 PetscCall(MatGetType(S_all, &Stype)); 1083 1084 /* we can reuse the solvers if we are not using the economic version */ 1085 reuse_solvers = (PetscBool)(reuse_solvers && !economic); 1086 if (!sub_schurs->gdsw) { 1087 factor_workaround = (PetscBool)(reuse_solvers && factor_workaround); 1088 if (!sub_schurs->is_posdef && factor_workaround && compute_Stilda && size_active_schur) reuse_solvers = factor_workaround = PETSC_FALSE; 1089 } 1090 solver_S = PETSC_TRUE; 1091 1092 /* update the Schur complement with the change of basis on the pressures */ 1093 if (benign_n) { 1094 const PetscScalar *cs_AIB; 1095 PetscScalar *S_data, *AIIm1_data; 1096 Mat S2 = NULL, S3 = NULL; /* dbg */ 1097 PetscScalar *S2_data, *S3_data; /* dbg */ 1098 Vec v, benign_AIIm1_ones; 1099 PetscInt sizeA; 1100 1101 PetscCall(MatDenseGetArray(S_all, &S_data)); 1102 PetscCall(MatCreateVecs(A, &v, &benign_AIIm1_ones)); 1103 PetscCall(VecGetSize(v, &sizeA)); 1104 PetscCall(MatMumpsSetIcntl(F, 26, 0)); 1105 #if defined(PETSC_HAVE_MKL_PARDISO) 1106 PetscCall(MatMkl_PardisoSetCntl(F, 70, 1)); 1107 #endif 1108 PetscCall(MatDenseGetArrayRead(cs_AIB_mat, &cs_AIB)); 1109 PetscCall(MatDenseGetArray(benign_AIIm1_ones_mat, &AIIm1_data)); 1110 if (matl_dbg_viewer) { 1111 PetscCall(MatDuplicate(S_all, MAT_DO_NOT_COPY_VALUES, &S2)); 1112 PetscCall(MatDuplicate(S_all, MAT_DO_NOT_COPY_VALUES, &S3)); 1113 PetscCall(MatDenseGetArray(S2, &S2_data)); 1114 PetscCall(MatDenseGetArray(S3, &S3_data)); 1115 } 1116 for (i = 0; i < benign_n; i++) { 1117 PetscScalar *array, sum = 0., one = 1., *sums; 1118 const PetscInt *idxs; 1119 PetscInt k, j, nz; 1120 PetscBLASInt B_k, B_n; 1121 1122 PetscCall(PetscCalloc1(benign_n, &sums)); 1123 PetscCall(VecPlaceArray(benign_AIIm1_ones, AIIm1_data + sizeA * i)); 1124 PetscCall(VecCopy(benign_AIIm1_ones, v)); 1125 PetscCall(MatSolve(F, v, benign_AIIm1_ones)); 1126 PetscCall(MatMult(A, benign_AIIm1_ones, v)); 1127 PetscCall(VecResetArray(benign_AIIm1_ones)); 1128 /* p0 dofs (eliminated) are excluded from the sums */ 1129 for (k = 0; k < benign_n; k++) { 1130 PetscCall(ISGetLocalSize(is_p_r[k], &nz)); 1131 PetscCall(ISGetIndices(is_p_r[k], &idxs)); 1132 for (j = 0; j < nz - 1; j++) sums[k] -= AIIm1_data[idxs[j] + sizeA * i]; 1133 PetscCall(ISRestoreIndices(is_p_r[k], &idxs)); 1134 } 1135 PetscCall(VecGetArrayRead(v, (const PetscScalar **)&array)); 1136 if (matl_dbg_viewer) { 1137 Vec vv; 1138 char name[16]; 1139 1140 PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, size_schur, array + n_I, &vv)); 1141 PetscCall(PetscSNPrintf(name, sizeof(name), "Pvs%" PetscInt_FMT, i)); 1142 PetscCall(PetscObjectSetName((PetscObject)vv, name)); 1143 PetscCall(VecView(vv, matl_dbg_viewer)); 1144 } 1145 /* perform sparse rank updates on symmetric Schur (TODO: move outside of the loop?) */ 1146 /* cs_AIB already scaled by 1./nz */ 1147 B_k = 1; 1148 PetscCall(PetscBLASIntCast(size_schur, &B_n)); 1149 for (k = 0; k < benign_n; k++) { 1150 sum = sums[k]; 1151 1152 if (PetscAbsScalar(sum) == 0.0) continue; 1153 if (k == i) { 1154 if (B_n) PetscCallBLAS("BLASsyrk", BLASsyrk_("L", "N", &B_n, &B_k, &sum, cs_AIB + i * size_schur, &B_n, &one, S_data, &B_n)); 1155 if (matl_dbg_viewer && B_n) PetscCallBLAS("BLASsyrk", BLASsyrk_("L", "N", &B_n, &B_k, &sum, cs_AIB + i * size_schur, &B_n, &one, S3_data, &B_n)); 1156 } else { /* XXX Is it correct to use symmetric rank-2 update with half of the sum? */ 1157 sum /= 2.0; 1158 if (B_n) PetscCallBLAS("BLASsyr2k", BLASsyr2k_("L", "N", &B_n, &B_k, &sum, cs_AIB + k * size_schur, &B_n, cs_AIB + i * size_schur, &B_n, &one, S_data, &B_n)); 1159 if (matl_dbg_viewer && B_n) PetscCallBLAS("BLASsyr2k", BLASsyr2k_("L", "N", &B_n, &B_k, &sum, cs_AIB + k * size_schur, &B_n, cs_AIB + i * size_schur, &B_n, &one, S3_data, &B_n)); 1160 } 1161 } 1162 sum = 1.; 1163 if (B_n) PetscCallBLAS("BLASsyr2k", BLASsyr2k_("L", "N", &B_n, &B_k, &sum, array + n_I, &B_n, cs_AIB + i * size_schur, &B_n, &one, S_data, &B_n)); 1164 if (matl_dbg_viewer && B_n) PetscCallBLAS("BLASsyr2k", BLASsyr2k_("L", "N", &B_n, &B_k, &sum, array + n_I, &B_n, cs_AIB + i * size_schur, &B_n, &one, S2_data, &B_n)); 1165 PetscCall(VecRestoreArrayRead(v, (const PetscScalar **)&array)); 1166 /* set p0 entry of AIIm1_ones to zero */ 1167 PetscCall(ISGetLocalSize(is_p_r[i], &nz)); 1168 PetscCall(ISGetIndices(is_p_r[i], &idxs)); 1169 for (j = 0; j < benign_n; j++) AIIm1_data[idxs[nz - 1] + sizeA * j] = 0.; 1170 PetscCall(ISRestoreIndices(is_p_r[i], &idxs)); 1171 PetscCall(PetscFree(sums)); 1172 } 1173 PetscCall(VecDestroy(&benign_AIIm1_ones)); 1174 if (matl_dbg_viewer) { 1175 PetscCall(MatDenseRestoreArray(S2, &S2_data)); 1176 PetscCall(MatDenseRestoreArray(S3, &S3_data)); 1177 } 1178 if (!S_lower_triangular) { /* I need to expand the upper triangular data (column-oriented) */ 1179 PetscInt k, j; 1180 for (k = 0; k < size_schur; k++) { 1181 for (j = k; j < size_schur; j++) S_data[j * size_schur + k] = PetscConj(S_data[k * size_schur + j]); 1182 } 1183 } 1184 1185 /* restore defaults */ 1186 PetscCall(MatMumpsSetIcntl(F, 26, -1)); 1187 #if defined(PETSC_HAVE_MKL_PARDISO) 1188 PetscCall(MatMkl_PardisoSetCntl(F, 70, 0)); 1189 #endif 1190 PetscCall(MatDenseRestoreArrayRead(cs_AIB_mat, &cs_AIB)); 1191 PetscCall(MatDenseRestoreArray(benign_AIIm1_ones_mat, &AIIm1_data)); 1192 PetscCall(VecDestroy(&v)); 1193 PetscCall(MatDenseRestoreArray(S_all, &S_data)); 1194 if (matl_dbg_viewer) { 1195 Mat S; 1196 1197 PetscCall(MatFactorRestoreSchurComplement(F, &S_all, MAT_FACTOR_SCHUR_UNFACTORED)); 1198 PetscCall(MatFactorCreateSchurComplement(F, &S, NULL)); 1199 PetscCall(PetscObjectSetName((PetscObject)S, "Sb")); 1200 PetscCall(MatView(S, matl_dbg_viewer)); 1201 PetscCall(MatDestroy(&S)); 1202 PetscCall(PetscObjectSetName((PetscObject)S2, "S2P")); 1203 PetscCall(MatView(S2, matl_dbg_viewer)); 1204 PetscCall(PetscObjectSetName((PetscObject)S3, "S3P")); 1205 PetscCall(MatView(S3, matl_dbg_viewer)); 1206 PetscCall(PetscObjectSetName((PetscObject)cs_AIB_mat, "cs")); 1207 PetscCall(MatView(cs_AIB_mat, matl_dbg_viewer)); 1208 PetscCall(MatFactorGetSchurComplement(F, &S_all, NULL)); 1209 } 1210 PetscCall(MatDestroy(&S2)); 1211 PetscCall(MatDestroy(&S3)); 1212 } 1213 if (!reuse_solvers) { 1214 for (i = 0; i < benign_n; i++) PetscCall(ISDestroy(&is_p_r[i])); 1215 PetscCall(PetscFree(is_p_r)); 1216 PetscCall(MatDestroy(&cs_AIB_mat)); 1217 PetscCall(MatDestroy(&benign_AIIm1_ones_mat)); 1218 } 1219 } else { /* we can't use MatFactor when size_schur == size_of_the_problem */ 1220 PetscCall(MatConvert(A, MATSEQDENSE, MAT_INITIAL_MATRIX, &S_all)); 1221 PetscCall(MatGetType(S_all, &Stype)); 1222 reuse_solvers = PETSC_FALSE; /* TODO: why we can't reuse the solvers here? */ 1223 factor_workaround = PETSC_FALSE; 1224 solver_S = PETSC_FALSE; 1225 } 1226 1227 if (reuse_solvers) { 1228 Mat A_II, pA_II, Afake; 1229 Vec vec1_B; 1230 PCBDDCReuseSolvers msolv_ctx; 1231 PetscInt n_R; 1232 1233 if (sub_schurs->reuse_solver) { 1234 PetscCall(PCBDDCReuseSolversReset(sub_schurs->reuse_solver)); 1235 } else { 1236 PetscCall(PetscNew(&sub_schurs->reuse_solver)); 1237 } 1238 msolv_ctx = sub_schurs->reuse_solver; 1239 PetscCall(MatSchurComplementGetSubMatrices(sub_schurs->S, &A_II, &pA_II, NULL, NULL, NULL)); 1240 PetscCall(PetscObjectReference((PetscObject)F)); 1241 msolv_ctx->F = F; 1242 PetscCall(MatCreateVecs(F, &msolv_ctx->sol, NULL)); 1243 /* currently PETSc has no support for MatSolve(F,x,x), so cheat and let rhs and sol share the same memory */ 1244 { 1245 PetscScalar *array; 1246 PetscInt n; 1247 1248 PetscCall(VecGetLocalSize(msolv_ctx->sol, &n)); 1249 PetscCall(VecGetArray(msolv_ctx->sol, &array)); 1250 PetscCall(VecCreateSeqWithArray(PetscObjectComm((PetscObject)msolv_ctx->sol), 1, n, array, &msolv_ctx->rhs)); 1251 PetscCall(VecRestoreArray(msolv_ctx->sol, &array)); 1252 } 1253 msolv_ctx->has_vertices = schur_has_vertices; 1254 1255 /* interior solver */ 1256 PetscCall(PCCreate(PetscObjectComm((PetscObject)A_II), &msolv_ctx->interior_solver)); 1257 PetscCall(PCSetOperators(msolv_ctx->interior_solver, A_II, pA_II)); 1258 PetscCall(PCSetType(msolv_ctx->interior_solver, PCSHELL)); 1259 PetscCall(PCShellSetName(msolv_ctx->interior_solver, "Interior solver (w/o Schur factorization)")); 1260 PetscCall(PCShellSetContext(msolv_ctx->interior_solver, msolv_ctx)); 1261 PetscCall(PCShellSetView(msolv_ctx->interior_solver, PCBDDCReuseSolvers_View)); 1262 PetscCall(PCShellSetApply(msolv_ctx->interior_solver, PCBDDCReuseSolvers_Interior)); 1263 PetscCall(PCShellSetApplyTranspose(msolv_ctx->interior_solver, PCBDDCReuseSolvers_InteriorTranspose)); 1264 if (sub_schurs->gdsw) PetscCall(PCShellSetDestroy(msolv_ctx->interior_solver, PCBDDCReuseSolvers_Destroy)); 1265 1266 /* correction solver */ 1267 if (!sub_schurs->gdsw) { 1268 PetscCall(PCCreate(PetscObjectComm((PetscObject)A_II), &msolv_ctx->correction_solver)); 1269 PetscCall(PCSetType(msolv_ctx->correction_solver, PCSHELL)); 1270 PetscCall(PCShellSetName(msolv_ctx->correction_solver, "Correction solver (with Schur factorization)")); 1271 PetscCall(PCShellSetContext(msolv_ctx->correction_solver, msolv_ctx)); 1272 PetscCall(PCShellSetView(msolv_ctx->interior_solver, PCBDDCReuseSolvers_View)); 1273 PetscCall(PCShellSetApply(msolv_ctx->correction_solver, PCBDDCReuseSolvers_Correction)); 1274 PetscCall(PCShellSetApplyTranspose(msolv_ctx->correction_solver, PCBDDCReuseSolvers_CorrectionTranspose)); 1275 1276 /* scatter and vecs for Schur complement solver */ 1277 PetscCall(MatCreateVecs(S_all, &msolv_ctx->sol_B, &msolv_ctx->rhs_B)); 1278 PetscCall(MatCreateVecs(sub_schurs->S, &vec1_B, NULL)); 1279 if (!schur_has_vertices) { 1280 PetscCall(ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap, IS_GTOLM_DROP, is_A_all, &msolv_ctx->is_B)); 1281 PetscCall(VecScatterCreate(vec1_B, msolv_ctx->is_B, msolv_ctx->sol_B, NULL, &msolv_ctx->correction_scatter_B)); 1282 PetscCall(PetscObjectReference((PetscObject)is_A_all)); 1283 msolv_ctx->is_R = is_A_all; 1284 } else { 1285 IS is_B_all; 1286 const PetscInt *idxs; 1287 PetscInt dual, n_v, n; 1288 1289 PetscCall(ISGetLocalSize(sub_schurs->is_vertices, &n_v)); 1290 dual = size_schur - n_v; 1291 PetscCall(ISGetLocalSize(is_A_all, &n)); 1292 PetscCall(ISGetIndices(is_A_all, &idxs)); 1293 PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)is_A_all), dual, idxs + n_I, PETSC_COPY_VALUES, &is_B_all)); 1294 PetscCall(ISGlobalToLocalMappingApplyIS(sub_schurs->BtoNmap, IS_GTOLM_DROP, is_B_all, &msolv_ctx->is_B)); 1295 PetscCall(ISDestroy(&is_B_all)); 1296 PetscCall(ISCreateStride(PetscObjectComm((PetscObject)is_A_all), dual, 0, 1, &is_B_all)); 1297 PetscCall(VecScatterCreate(vec1_B, msolv_ctx->is_B, msolv_ctx->sol_B, is_B_all, &msolv_ctx->correction_scatter_B)); 1298 PetscCall(ISDestroy(&is_B_all)); 1299 PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)is_A_all), n - n_v, idxs, PETSC_COPY_VALUES, &msolv_ctx->is_R)); 1300 PetscCall(ISRestoreIndices(is_A_all, &idxs)); 1301 } 1302 PetscCall(ISGetLocalSize(msolv_ctx->is_R, &n_R)); 1303 PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF, n_R, n_R, 0, NULL, &Afake)); 1304 PetscCall(MatAssemblyBegin(Afake, MAT_FINAL_ASSEMBLY)); 1305 PetscCall(MatAssemblyEnd(Afake, MAT_FINAL_ASSEMBLY)); 1306 PetscCall(PCSetOperators(msolv_ctx->correction_solver, Afake, Afake)); 1307 PetscCall(MatDestroy(&Afake)); 1308 PetscCall(VecDestroy(&vec1_B)); 1309 } 1310 /* communicate benign info to solver context */ 1311 if (benign_n) { 1312 PetscScalar *array; 1313 1314 msolv_ctx->benign_n = benign_n; 1315 msolv_ctx->benign_zerodiag_subs = is_p_r; 1316 PetscCall(PetscMalloc1(benign_n, &msolv_ctx->benign_save_vals)); 1317 msolv_ctx->benign_csAIB = cs_AIB_mat; 1318 PetscCall(MatCreateVecs(cs_AIB_mat, &msolv_ctx->benign_corr_work, NULL)); 1319 PetscCall(VecGetArray(msolv_ctx->benign_corr_work, &array)); 1320 PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, size_schur, array, &msolv_ctx->benign_dummy_schur_vec)); 1321 PetscCall(VecRestoreArray(msolv_ctx->benign_corr_work, &array)); 1322 msolv_ctx->benign_AIIm1ones = benign_AIIm1_ones_mat; 1323 } 1324 } else { 1325 if (sub_schurs->reuse_solver) PetscCall(PCBDDCReuseSolversReset(sub_schurs->reuse_solver)); 1326 PetscCall(PetscFree(sub_schurs->reuse_solver)); 1327 } 1328 PetscCall(MatDestroy(&A)); 1329 PetscCall(ISDestroy(&is_A_all)); 1330 1331 /* Work arrays */ 1332 PetscCall(PetscMalloc1(max_subset_size * max_subset_size, &work)); 1333 1334 /* S_Ej_all */ 1335 cum = cum2 = 0; 1336 PetscCall(MatDenseGetArrayRead(S_all, &rS_data)); 1337 PetscCall(MatSeqAIJGetArray(sub_schurs->S_Ej_all, &SEj_arr)); 1338 if (sub_schurs->sum_S_Ej_inv_all) PetscCall(MatSeqAIJGetArray(sub_schurs->sum_S_Ej_inv_all, &SEjinv_arr)); 1339 if (sub_schurs->gdsw) PetscCall(MatCreateSubMatrices(sub_schurs->A, sub_schurs->n_subs, sub_schurs->is_subs, sub_schurs->is_subs, MAT_INITIAL_MATRIX, &gdswA)); 1340 for (i = 0; i < sub_schurs->n_subs; i++) { 1341 PetscInt j; 1342 1343 /* get S_E (or K^i_EE for GDSW) */ 1344 PetscCall(ISGetLocalSize(sub_schurs->is_subs[i], &subset_size)); 1345 if (sub_schurs->gdsw) { 1346 Mat T; 1347 1348 PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, subset_size, subset_size, work, &T)); 1349 PetscCall(MatConvert(gdswA[i], MATDENSE, MAT_REUSE_MATRIX, &T)); 1350 PetscCall(MatDestroy(&T)); 1351 } else { 1352 if (S_lower_triangular) { /* I need to expand the upper triangular data (column-oriented) */ 1353 PetscInt k; 1354 for (k = 0; k < subset_size; k++) { 1355 for (j = k; j < subset_size; j++) { 1356 work[k * subset_size + j] = rS_data[cum2 + k * size_schur + j]; 1357 work[j * subset_size + k] = PetscConj(rS_data[cum2 + k * size_schur + j]); 1358 } 1359 } 1360 } else { /* just copy to workspace */ 1361 PetscInt k; 1362 for (k = 0; k < subset_size; k++) { 1363 for (j = 0; j < subset_size; j++) work[k * subset_size + j] = rS_data[cum2 + k * size_schur + j]; 1364 } 1365 } 1366 } 1367 /* insert S_E values */ 1368 if (sub_schurs->change) { 1369 Mat change_sub, SEj, T; 1370 1371 /* change basis */ 1372 PetscCall(KSPGetOperators(sub_schurs->change[i], &change_sub, NULL)); 1373 PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, subset_size, subset_size, work, &SEj)); 1374 if (!sub_schurs->change_with_qr) { /* currently there's no support for PtAP with P SeqAIJ */ 1375 Mat T2; 1376 PetscCall(MatTransposeMatMult(change_sub, SEj, MAT_INITIAL_MATRIX, 1.0, &T2)); 1377 PetscCall(MatMatMult(T2, change_sub, MAT_INITIAL_MATRIX, 1.0, &T)); 1378 PetscCall(MatConvert(T, MATSEQDENSE, MAT_INPLACE_MATRIX, &T)); 1379 PetscCall(MatDestroy(&T2)); 1380 } else { 1381 PetscCall(MatPtAP(SEj, change_sub, MAT_INITIAL_MATRIX, 1.0, &T)); 1382 } 1383 PetscCall(MatCopy(T, SEj, SAME_NONZERO_PATTERN)); 1384 PetscCall(MatDestroy(&T)); 1385 PetscCall(MatZeroRowsColumnsIS(SEj, sub_schurs->change_primal_sub[i], 1.0, NULL, NULL)); 1386 PetscCall(MatDestroy(&SEj)); 1387 } 1388 PetscCall(PetscArraycpy(SEj_arr, work, subset_size * subset_size)); 1389 if (compute_Stilda) { 1390 if (deluxe) { /* if adaptivity is requested, invert S_E blocks */ 1391 Mat M; 1392 const PetscScalar *vals; 1393 PetscBool isdense, isdensecuda; 1394 1395 PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, subset_size, subset_size, work, &M)); 1396 PetscCall(MatSetOption(M, MAT_SPD, sub_schurs->is_posdef)); 1397 PetscCall(MatSetOption(M, MAT_HERMITIAN, sub_schurs->is_hermitian)); 1398 if (!PetscBTLookup(sub_schurs->is_edge, i)) PetscCall(MatSetType(M, Stype)); 1399 PetscCall(PetscObjectTypeCompare((PetscObject)M, MATSEQDENSE, &isdense)); 1400 PetscCall(PetscObjectTypeCompare((PetscObject)M, MATSEQDENSECUDA, &isdensecuda)); 1401 switch (sub_schurs->mat_factor_type) { 1402 case MAT_FACTOR_CHOLESKY: 1403 PetscCall(MatCholeskyFactor(M, NULL, NULL)); 1404 break; 1405 case MAT_FACTOR_LU: 1406 PetscCall(MatLUFactor(M, NULL, NULL, NULL)); 1407 break; 1408 default: 1409 SETERRQ(PetscObjectComm((PetscObject)F), PETSC_ERR_SUP, "Unsupported factor type %s", MatFactorTypes[sub_schurs->mat_factor_type]); 1410 } 1411 if (isdense) { 1412 PetscCall(MatSeqDenseInvertFactors_Private(M)); 1413 #if defined(PETSC_HAVE_CUDA) 1414 } else if (isdensecuda) { 1415 PetscCall(MatSeqDenseCUDAInvertFactors_Internal(M)); 1416 #endif 1417 } else SETERRQ(PetscObjectComm((PetscObject)M), PETSC_ERR_SUP, "Not implemented for type %s", Stype); 1418 PetscCall(MatDenseGetArrayRead(M, &vals)); 1419 PetscCall(PetscArraycpy(SEjinv_arr, vals, subset_size * subset_size)); 1420 PetscCall(MatDenseRestoreArrayRead(M, &vals)); 1421 PetscCall(MatDestroy(&M)); 1422 } else if (scaling) { /* not using deluxe */ 1423 Mat SEj; 1424 Vec D; 1425 PetscScalar *array; 1426 1427 PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, subset_size, subset_size, work, &SEj)); 1428 PetscCall(VecGetArray(Dall, &array)); 1429 PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, subset_size, array + cum, &D)); 1430 PetscCall(VecRestoreArray(Dall, &array)); 1431 PetscCall(VecShift(D, -1.)); 1432 PetscCall(MatDiagonalScale(SEj, D, D)); 1433 PetscCall(MatDestroy(&SEj)); 1434 PetscCall(VecDestroy(&D)); 1435 PetscCall(PetscArraycpy(SEj_arr, work, subset_size * subset_size)); 1436 } 1437 } 1438 cum += subset_size; 1439 cum2 += subset_size * (size_schur + 1); 1440 SEj_arr += subset_size * subset_size; 1441 if (SEjinv_arr) SEjinv_arr += subset_size * subset_size; 1442 } 1443 if (sub_schurs->gdsw) PetscCall(MatDestroySubMatrices(sub_schurs->n_subs, &gdswA)); 1444 PetscCall(MatDenseRestoreArrayRead(S_all, &rS_data)); 1445 PetscCall(MatSeqAIJRestoreArray(sub_schurs->S_Ej_all, &SEj_arr)); 1446 if (sub_schurs->sum_S_Ej_inv_all) PetscCall(MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_inv_all, &SEjinv_arr)); 1447 if (solver_S) PetscCall(MatFactorRestoreSchurComplement(F, &S_all, MAT_FACTOR_SCHUR_UNFACTORED)); 1448 1449 /* may prevent from unneeded copies, since MUMPS or MKL_Pardiso always use CPU memory 1450 however, preliminary tests indicate using GPUs is still faster in the solve phase */ 1451 #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA) 1452 if (reuse_solvers) { 1453 Mat St; 1454 MatFactorSchurStatus st; 1455 1456 flg = PETSC_FALSE; 1457 PetscCall(PetscOptionsGetBool(NULL, sub_schurs->prefix, "-sub_schurs_schur_pin_to_cpu", &flg, NULL)); 1458 PetscCall(MatFactorGetSchurComplement(F, &St, &st)); 1459 PetscCall(MatBindToCPU(St, flg)); 1460 PetscCall(MatFactorRestoreSchurComplement(F, &St, st)); 1461 } 1462 #endif 1463 1464 schur_factor = NULL; 1465 if (compute_Stilda && size_active_schur) { 1466 if (sub_schurs->n_subs == 1 && size_schur == size_active_schur && deluxe) { /* we already computed the inverse */ 1467 PetscCall(MatSeqAIJGetArrayWrite(sub_schurs->sum_S_Ej_tilda_all, &SEjinv_arr)); 1468 PetscCall(PetscArraycpy(SEjinv_arr, work, size_schur * size_schur)); 1469 PetscCall(MatSeqAIJRestoreArrayWrite(sub_schurs->sum_S_Ej_tilda_all, &SEjinv_arr)); 1470 } else { 1471 Mat S_all_inv = NULL; 1472 1473 if (solver_S && !sub_schurs->gdsw) { 1474 /* for adaptive selection we need S^-1; for solver reusage we need S_\Delta\Delta^-1. 1475 The latter is not the principal subminor for S^-1. However, the factors can be reused since S_\Delta\Delta is the leading principal submatrix of S */ 1476 if (factor_workaround) { /* invert without calling MatFactorInvertSchurComplement, since we are hacking */ 1477 PetscScalar *data; 1478 PetscInt nd = 0; 1479 1480 PetscCheck(use_potr, PETSC_COMM_SELF, PETSC_ERR_SUP, "Factor update not yet implemented for non SPD matrices"); 1481 PetscCall(MatFactorGetSchurComplement(F, &S_all_inv, NULL)); 1482 PetscCall(MatDenseGetArray(S_all_inv, &data)); 1483 if (sub_schurs->is_dir) { /* dirichlet dofs could have different scalings */ 1484 PetscCall(ISGetLocalSize(sub_schurs->is_dir, &nd)); 1485 } 1486 1487 /* factor and invert activedofs and vertices (dirichlet dofs does not contribute) */ 1488 if (schur_has_vertices) { 1489 Mat M; 1490 PetscScalar *tdata; 1491 PetscInt nv = 0, news; 1492 1493 PetscCall(ISGetLocalSize(sub_schurs->is_vertices, &nv)); 1494 news = size_active_schur + nv; 1495 PetscCall(PetscCalloc1(news * news, &tdata)); 1496 for (i = 0; i < size_active_schur; i++) { 1497 PetscCall(PetscArraycpy(tdata + i * (news + 1), data + i * (size_schur + 1), size_active_schur - i)); 1498 PetscCall(PetscArraycpy(tdata + i * (news + 1) + size_active_schur - i, data + i * size_schur + size_active_schur + nd, nv)); 1499 } 1500 for (i = 0; i < nv; i++) { 1501 PetscInt k = i + size_active_schur; 1502 PetscCall(PetscArraycpy(tdata + k * (news + 1), data + (k + nd) * (size_schur + 1), nv - i)); 1503 } 1504 1505 PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, news, news, tdata, &M)); 1506 PetscCall(MatSetOption(M, MAT_SPD, PETSC_TRUE)); 1507 PetscCall(MatCholeskyFactor(M, NULL, NULL)); 1508 /* save the factors */ 1509 cum = 0; 1510 PetscCall(PetscMalloc1((size_active_schur * (size_active_schur + 1)) / 2 + nd, &schur_factor)); 1511 for (i = 0; i < size_active_schur; i++) { 1512 PetscCall(PetscArraycpy(schur_factor + cum, tdata + i * (news + 1), size_active_schur - i)); 1513 cum += size_active_schur - i; 1514 } 1515 for (i = 0; i < nd; i++) schur_factor[cum + i] = PetscSqrtReal(PetscRealPart(data[(i + size_active_schur) * (size_schur + 1)])); 1516 PetscCall(MatSeqDenseInvertFactors_Private(M)); 1517 /* move back just the active dofs to the Schur complement */ 1518 for (i = 0; i < size_active_schur; i++) PetscCall(PetscArraycpy(data + i * size_schur, tdata + i * news, size_active_schur)); 1519 PetscCall(PetscFree(tdata)); 1520 PetscCall(MatDestroy(&M)); 1521 } else { /* we can factorize and invert just the activedofs */ 1522 Mat M; 1523 PetscScalar *aux; 1524 1525 PetscCall(PetscMalloc1(nd, &aux)); 1526 for (i = 0; i < nd; i++) aux[i] = 1.0 / data[(i + size_active_schur) * (size_schur + 1)]; 1527 PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, size_active_schur, size_active_schur, data, &M)); 1528 PetscCall(MatDenseSetLDA(M, size_schur)); 1529 PetscCall(MatSetOption(M, MAT_SPD, PETSC_TRUE)); 1530 PetscCall(MatCholeskyFactor(M, NULL, NULL)); 1531 PetscCall(MatSeqDenseInvertFactors_Private(M)); 1532 PetscCall(MatDestroy(&M)); 1533 PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, size_schur, nd, data + size_active_schur * size_schur, &M)); 1534 PetscCall(MatZeroEntries(M)); 1535 PetscCall(MatDestroy(&M)); 1536 PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, nd, size_schur, data + size_active_schur, &M)); 1537 PetscCall(MatDenseSetLDA(M, size_schur)); 1538 PetscCall(MatZeroEntries(M)); 1539 PetscCall(MatDestroy(&M)); 1540 for (i = 0; i < nd; i++) data[(i + size_active_schur) * (size_schur + 1)] = aux[i]; 1541 PetscCall(PetscFree(aux)); 1542 } 1543 PetscCall(MatDenseRestoreArray(S_all_inv, &data)); 1544 } else { /* use MatFactor calls to invert S */ 1545 PetscCall(MatFactorInvertSchurComplement(F)); 1546 PetscCall(MatFactorGetSchurComplement(F, &S_all_inv, NULL)); 1547 } 1548 } else if (!sub_schurs->gdsw) { /* we need to invert explicitly since we are not using MatFactor for S */ 1549 PetscCall(PetscObjectReference((PetscObject)S_all)); 1550 S_all_inv = S_all; 1551 PetscCall(MatDenseGetArray(S_all_inv, &S_data)); 1552 PetscCall(PetscBLASIntCast(size_schur, &B_N)); 1553 PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); 1554 if (use_potr) { 1555 PetscCallBLAS("LAPACKpotrf", LAPACKpotrf_("L", &B_N, S_data, &B_N, &B_ierr)); 1556 PetscCheck(!B_ierr, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in POTRF Lapack routine %" PetscBLASInt_FMT, B_ierr); 1557 PetscCallBLAS("LAPACKpotri", LAPACKpotri_("L", &B_N, S_data, &B_N, &B_ierr)); 1558 PetscCheck(!B_ierr, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in POTRI Lapack routine %" PetscBLASInt_FMT, B_ierr); 1559 } else if (use_sytr) { 1560 PetscCallBLAS("LAPACKsytrf", LAPACKsytrf_("L", &B_N, S_data, &B_N, pivots, Bwork, &B_lwork, &B_ierr)); 1561 PetscCheck(!B_ierr, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in SYTRF Lapack routine %" PetscBLASInt_FMT, B_ierr); 1562 PetscCallBLAS("LAPACKsytri", LAPACKsytri_("L", &B_N, S_data, &B_N, pivots, Bwork, &B_ierr)); 1563 PetscCheck(!B_ierr, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in SYTRI Lapack routine %" PetscBLASInt_FMT, B_ierr); 1564 } else { 1565 PetscCallBLAS("LAPACKgetrf", LAPACKgetrf_(&B_N, &B_N, S_data, &B_N, pivots, &B_ierr)); 1566 PetscCheck(!B_ierr, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in GETRF Lapack routine %" PetscBLASInt_FMT, B_ierr); 1567 PetscCallBLAS("LAPACKgetri", LAPACKgetri_(&B_N, S_data, &B_N, pivots, Bwork, &B_lwork, &B_ierr)); 1568 PetscCheck(!B_ierr, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in GETRI Lapack routine %" PetscBLASInt_FMT, B_ierr); 1569 } 1570 PetscCall(PetscLogFlops(1.0 * size_schur * size_schur * size_schur)); 1571 PetscCall(PetscFPTrapPop()); 1572 PetscCall(MatDenseRestoreArray(S_all_inv, &S_data)); 1573 } else if (sub_schurs->gdsw) { 1574 Mat tS, tX, SEj, S_II, S_IE, S_EE; 1575 KSP pS_II; 1576 PC pS_II_pc; 1577 IS EE, II; 1578 PetscInt nS; 1579 1580 PetscCall(MatFactorCreateSchurComplement(F, &tS, NULL)); 1581 PetscCall(MatGetSize(tS, &nS, NULL)); 1582 PetscCall(MatSeqAIJGetArray(sub_schurs->sum_S_Ej_tilda_all, &SEjinv_arr)); 1583 for (i = 0, cum = 0; i < sub_schurs->n_subs; i++) { /* naive implementation */ 1584 PetscCall(ISGetLocalSize(sub_schurs->is_subs[i], &subset_size)); 1585 PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, subset_size, subset_size, SEjinv_arr, &SEj)); 1586 1587 PetscCall(ISCreateStride(PETSC_COMM_SELF, subset_size, cum, 1, &EE)); 1588 PetscCall(ISComplement(EE, 0, nS, &II)); 1589 PetscCall(MatCreateSubMatrix(tS, II, II, MAT_INITIAL_MATRIX, &S_II)); 1590 PetscCall(MatCreateSubMatrix(tS, II, EE, MAT_INITIAL_MATRIX, &S_IE)); 1591 PetscCall(MatCreateSubMatrix(tS, EE, EE, MAT_INITIAL_MATRIX, &S_EE)); 1592 PetscCall(ISDestroy(&II)); 1593 PetscCall(ISDestroy(&EE)); 1594 1595 PetscCall(KSPCreate(PETSC_COMM_SELF, &pS_II)); 1596 PetscCall(KSPSetNestLevel(pS_II, 1)); /* do not have direct access to a PC to provide the level of nesting of the KSP */ 1597 PetscCall(KSPSetType(pS_II, KSPPREONLY)); 1598 PetscCall(KSPGetPC(pS_II, &pS_II_pc)); 1599 PetscCall(PCSetType(pS_II_pc, PCSVD)); 1600 PetscCall(KSPSetOptionsPrefix(pS_II, sub_schurs->prefix)); 1601 PetscCall(KSPAppendOptionsPrefix(pS_II, "pseudo_")); 1602 PetscCall(KSPSetOperators(pS_II, S_II, S_II)); 1603 PetscCall(MatDestroy(&S_II)); 1604 PetscCall(KSPSetFromOptions(pS_II)); 1605 PetscCall(KSPSetUp(pS_II)); 1606 PetscCall(MatDuplicate(S_IE, MAT_DO_NOT_COPY_VALUES, &tX)); 1607 PetscCall(KSPMatSolve(pS_II, S_IE, tX)); 1608 PetscCall(KSPDestroy(&pS_II)); 1609 1610 PetscCall(MatTransposeMatMult(S_IE, tX, MAT_REUSE_MATRIX, PETSC_DETERMINE, &SEj)); 1611 PetscCall(MatDestroy(&S_IE)); 1612 PetscCall(MatDestroy(&tX)); 1613 PetscCall(MatAYPX(SEj, -1, S_EE, SAME_NONZERO_PATTERN)); 1614 PetscCall(MatDestroy(&S_EE)); 1615 1616 PetscCall(MatDestroy(&SEj)); 1617 cum += subset_size; 1618 SEjinv_arr += subset_size * subset_size; 1619 } 1620 PetscCall(MatDestroy(&tS)); 1621 PetscCall(MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_tilda_all, &SEjinv_arr)); 1622 } 1623 /* S_Ej_tilda_all */ 1624 cum = cum2 = 0; 1625 rS_data = NULL; 1626 if (S_all_inv) PetscCall(MatDenseGetArrayRead(S_all_inv, &rS_data)); 1627 PetscCall(MatSeqAIJGetArrayWrite(sub_schurs->sum_S_Ej_tilda_all, &SEjinv_arr)); 1628 for (i = 0; i < sub_schurs->n_subs; i++) { 1629 PetscInt j; 1630 1631 PetscCall(ISGetLocalSize(sub_schurs->is_subs[i], &subset_size)); 1632 /* get (St^-1)_E */ 1633 /* Unless we are changing the variables, I don't need to expand to upper triangular since St^-1 1634 will be properly accessed later during adaptive selection */ 1635 if (rS_data) { 1636 if (S_lower_triangular) { 1637 PetscInt k; 1638 if (sub_schurs->change) { 1639 for (k = 0; k < subset_size; k++) { 1640 for (j = k; j < subset_size; j++) { 1641 work[k * subset_size + j] = rS_data[cum2 + k * size_schur + j]; 1642 work[j * subset_size + k] = work[k * subset_size + j]; 1643 } 1644 } 1645 } else { 1646 for (k = 0; k < subset_size; k++) { 1647 for (j = k; j < subset_size; j++) work[k * subset_size + j] = rS_data[cum2 + k * size_schur + j]; 1648 } 1649 } 1650 } else { 1651 PetscInt k; 1652 for (k = 0; k < subset_size; k++) { 1653 for (j = 0; j < subset_size; j++) work[k * subset_size + j] = rS_data[cum2 + k * size_schur + j]; 1654 } 1655 } 1656 } 1657 if (sub_schurs->change) { 1658 Mat change_sub, SEj, T; 1659 PetscScalar val = sub_schurs->gdsw ? PETSC_SMALL : 1. / PETSC_SMALL; 1660 1661 /* change basis */ 1662 PetscCall(KSPGetOperators(sub_schurs->change[i], &change_sub, NULL)); 1663 PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, subset_size, subset_size, rS_data ? work : SEjinv_arr, &SEj)); 1664 if (!sub_schurs->change_with_qr) { /* currently there's no support for PtAP with P SeqAIJ */ 1665 Mat T2; 1666 PetscCall(MatTransposeMatMult(change_sub, SEj, MAT_INITIAL_MATRIX, 1.0, &T2)); 1667 PetscCall(MatMatMult(T2, change_sub, MAT_INITIAL_MATRIX, 1.0, &T)); 1668 PetscCall(MatDestroy(&T2)); 1669 PetscCall(MatConvert(T, MATSEQDENSE, MAT_INPLACE_MATRIX, &T)); 1670 } else { 1671 PetscCall(MatPtAP(SEj, change_sub, MAT_INITIAL_MATRIX, 1.0, &T)); 1672 } 1673 PetscCall(MatCopy(T, SEj, SAME_NONZERO_PATTERN)); 1674 PetscCall(MatDestroy(&T)); 1675 PetscCall(MatZeroRowsColumnsIS(SEj, sub_schurs->change_primal_sub[i], val, NULL, NULL)); 1676 PetscCall(MatDestroy(&SEj)); 1677 } 1678 if (rS_data) PetscCall(PetscArraycpy(SEjinv_arr, work, subset_size * subset_size)); 1679 cum += subset_size; 1680 cum2 += subset_size * (size_schur + 1); 1681 SEjinv_arr += subset_size * subset_size; 1682 } 1683 PetscCall(MatSeqAIJRestoreArrayWrite(sub_schurs->sum_S_Ej_tilda_all, &SEjinv_arr)); 1684 if (S_all_inv) { 1685 PetscCall(MatDenseRestoreArrayRead(S_all_inv, &rS_data)); 1686 if (solver_S) { 1687 if (schur_has_vertices) { 1688 PetscCall(MatFactorRestoreSchurComplement(F, &S_all_inv, MAT_FACTOR_SCHUR_FACTORED)); 1689 } else { 1690 PetscCall(MatFactorRestoreSchurComplement(F, &S_all_inv, MAT_FACTOR_SCHUR_INVERTED)); 1691 } 1692 } 1693 } 1694 PetscCall(MatDestroy(&S_all_inv)); 1695 } 1696 1697 /* move back factors if needed */ 1698 if (schur_has_vertices && factor_workaround && !sub_schurs->gdsw) { 1699 Mat S_tmp; 1700 PetscInt nd = 0; 1701 PetscScalar *data; 1702 1703 PetscCheck(use_potr, PETSC_COMM_SELF, PETSC_ERR_SUP, "Factor update not yet implemented for non SPD matrices"); 1704 PetscCheck(solver_S, PETSC_COMM_SELF, PETSC_ERR_PLIB, "This should not happen"); 1705 PetscCall(MatFactorGetSchurComplement(F, &S_tmp, NULL)); 1706 PetscCall(MatDenseGetArray(S_tmp, &data)); 1707 PetscCall(PetscArrayzero(data, size_schur * size_schur)); 1708 1709 if (S_lower_triangular) { 1710 cum = 0; 1711 for (i = 0; i < size_active_schur; i++) { 1712 PetscCall(PetscArraycpy(data + i * (size_schur + 1), schur_factor + cum, size_active_schur - i)); 1713 cum += size_active_schur - i; 1714 } 1715 } else { 1716 PetscCall(PetscArraycpy(data, schur_factor, size_schur * size_schur)); 1717 } 1718 if (sub_schurs->is_dir) { 1719 PetscCall(ISGetLocalSize(sub_schurs->is_dir, &nd)); 1720 for (i = 0; i < nd; i++) data[(i + size_active_schur) * (size_schur + 1)] = schur_factor[cum + i]; 1721 } 1722 /* workaround: since I cannot modify the matrices used inside the solvers for the forward and backward substitutions, 1723 set the diagonal entry of the Schur factor to a very large value */ 1724 for (i = size_active_schur + nd; i < size_schur; i++) data[i * (size_schur + 1)] = infty; 1725 PetscCall(MatDenseRestoreArray(S_tmp, &data)); 1726 PetscCall(MatFactorRestoreSchurComplement(F, &S_tmp, MAT_FACTOR_SCHUR_FACTORED)); 1727 } 1728 } else if (factor_workaround && !sub_schurs->gdsw) { /* we need to eliminate any unneeded coupling */ 1729 PetscScalar *data; 1730 PetscInt nd = 0; 1731 1732 if (sub_schurs->is_dir) { /* dirichlet dofs could have different scalings */ 1733 PetscCall(ISGetLocalSize(sub_schurs->is_dir, &nd)); 1734 } 1735 PetscCall(MatFactorGetSchurComplement(F, &S_all, NULL)); 1736 PetscCall(MatDenseGetArray(S_all, &data)); 1737 for (i = 0; i < size_active_schur; i++) PetscCall(PetscArrayzero(data + i * size_schur + size_active_schur, size_schur - size_active_schur)); 1738 for (i = size_active_schur + nd; i < size_schur; i++) { 1739 PetscCall(PetscArrayzero(data + i * size_schur + size_active_schur, size_schur - size_active_schur)); 1740 data[i * (size_schur + 1)] = infty; 1741 } 1742 PetscCall(MatDenseRestoreArray(S_all, &data)); 1743 PetscCall(MatFactorRestoreSchurComplement(F, &S_all, MAT_FACTOR_SCHUR_UNFACTORED)); 1744 } 1745 PetscCall(PetscFree(work)); 1746 PetscCall(PetscFree(schur_factor)); 1747 PetscCall(VecDestroy(&Dall)); 1748 } 1749 PetscCall(ISDestroy(&is_I_layer)); 1750 PetscCall(MatDestroy(&S_all)); 1751 PetscCall(MatDestroy(&A_BB)); 1752 PetscCall(MatDestroy(&A_IB)); 1753 PetscCall(MatDestroy(&A_BI)); 1754 PetscCall(MatDestroy(&F)); 1755 1756 PetscCall(PetscMalloc1(sub_schurs->n_subs, &nnz)); 1757 for (i = 0; i < sub_schurs->n_subs; i++) PetscCall(ISGetLocalSize(sub_schurs->is_subs[i], &nnz[i])); 1758 PetscCall(ISCreateGeneral(PETSC_COMM_SELF, sub_schurs->n_subs, nnz, PETSC_OWN_POINTER, &is_I_layer)); 1759 PetscCall(MatSetVariableBlockSizes(sub_schurs->S_Ej_all, sub_schurs->n_subs, nnz)); 1760 PetscCall(MatAssemblyBegin(sub_schurs->S_Ej_all, MAT_FINAL_ASSEMBLY)); 1761 PetscCall(MatAssemblyEnd(sub_schurs->S_Ej_all, MAT_FINAL_ASSEMBLY)); 1762 if (compute_Stilda) { 1763 PetscCall(MatSetVariableBlockSizes(sub_schurs->sum_S_Ej_tilda_all, sub_schurs->n_subs, nnz)); 1764 PetscCall(MatAssemblyBegin(sub_schurs->sum_S_Ej_tilda_all, MAT_FINAL_ASSEMBLY)); 1765 PetscCall(MatAssemblyEnd(sub_schurs->sum_S_Ej_tilda_all, MAT_FINAL_ASSEMBLY)); 1766 if (deluxe) { 1767 PetscCall(MatSetVariableBlockSizes(sub_schurs->sum_S_Ej_inv_all, sub_schurs->n_subs, nnz)); 1768 PetscCall(MatAssemblyBegin(sub_schurs->sum_S_Ej_inv_all, MAT_FINAL_ASSEMBLY)); 1769 PetscCall(MatAssemblyEnd(sub_schurs->sum_S_Ej_inv_all, MAT_FINAL_ASSEMBLY)); 1770 } 1771 } 1772 PetscCall(ISDestroy(&is_I_layer)); 1773 1774 /* Get local part of (\sum_j S_Ej) */ 1775 if (!sub_schurs->sum_S_Ej_all) PetscCall(MatDuplicate(sub_schurs->S_Ej_all, MAT_DO_NOT_COPY_VALUES, &sub_schurs->sum_S_Ej_all)); 1776 PetscCall(VecSet(gstash, 0.0)); 1777 PetscCall(MatSeqAIJGetArray(sub_schurs->S_Ej_all, &stasharray)); 1778 PetscCall(VecPlaceArray(lstash, stasharray)); 1779 PetscCall(VecScatterBegin(sstash, lstash, gstash, ADD_VALUES, SCATTER_FORWARD)); 1780 PetscCall(VecScatterEnd(sstash, lstash, gstash, ADD_VALUES, SCATTER_FORWARD)); 1781 PetscCall(MatSeqAIJRestoreArray(sub_schurs->S_Ej_all, &stasharray)); 1782 PetscCall(VecResetArray(lstash)); 1783 PetscCall(MatSeqAIJGetArray(sub_schurs->sum_S_Ej_all, &stasharray)); 1784 PetscCall(VecPlaceArray(lstash, stasharray)); 1785 PetscCall(VecScatterBegin(sstash, gstash, lstash, INSERT_VALUES, SCATTER_REVERSE)); 1786 PetscCall(VecScatterEnd(sstash, gstash, lstash, INSERT_VALUES, SCATTER_REVERSE)); 1787 PetscCall(MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_all, &stasharray)); 1788 PetscCall(VecResetArray(lstash)); 1789 1790 /* Get local part of (\sum_j S^-1_Ej) (\sum_j St^-1_Ej) */ 1791 if (compute_Stilda) { 1792 PetscCall(VecSet(gstash, 0.0)); 1793 PetscCall(MatSeqAIJGetArray(sub_schurs->sum_S_Ej_tilda_all, &stasharray)); 1794 PetscCall(VecPlaceArray(lstash, stasharray)); 1795 PetscCall(VecScatterBegin(sstash, lstash, gstash, ADD_VALUES, SCATTER_FORWARD)); 1796 PetscCall(VecScatterEnd(sstash, lstash, gstash, ADD_VALUES, SCATTER_FORWARD)); 1797 PetscCall(VecScatterBegin(sstash, gstash, lstash, INSERT_VALUES, SCATTER_REVERSE)); 1798 PetscCall(VecScatterEnd(sstash, gstash, lstash, INSERT_VALUES, SCATTER_REVERSE)); 1799 PetscCall(MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_tilda_all, &stasharray)); 1800 PetscCall(VecResetArray(lstash)); 1801 if (deluxe) { 1802 PetscCall(VecSet(gstash, 0.0)); 1803 PetscCall(MatSeqAIJGetArray(sub_schurs->sum_S_Ej_inv_all, &stasharray)); 1804 PetscCall(VecPlaceArray(lstash, stasharray)); 1805 PetscCall(VecScatterBegin(sstash, lstash, gstash, ADD_VALUES, SCATTER_FORWARD)); 1806 PetscCall(VecScatterEnd(sstash, lstash, gstash, ADD_VALUES, SCATTER_FORWARD)); 1807 PetscCall(VecScatterBegin(sstash, gstash, lstash, INSERT_VALUES, SCATTER_REVERSE)); 1808 PetscCall(VecScatterEnd(sstash, gstash, lstash, INSERT_VALUES, SCATTER_REVERSE)); 1809 PetscCall(MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_inv_all, &stasharray)); 1810 PetscCall(VecResetArray(lstash)); 1811 } else if (!sub_schurs->gdsw) { 1812 PetscScalar *array; 1813 PetscInt cum; 1814 1815 PetscCall(MatSeqAIJGetArray(sub_schurs->sum_S_Ej_tilda_all, &array)); 1816 cum = 0; 1817 for (i = 0; i < sub_schurs->n_subs; i++) { 1818 PetscCall(ISGetLocalSize(sub_schurs->is_subs[i], &subset_size)); 1819 PetscCall(PetscBLASIntCast(subset_size, &B_N)); 1820 PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); 1821 if (use_potr) { 1822 PetscCallBLAS("LAPACKpotrf", LAPACKpotrf_("L", &B_N, array + cum, &B_N, &B_ierr)); 1823 PetscCheck(!B_ierr, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in POTRF Lapack routine %" PetscBLASInt_FMT, B_ierr); 1824 PetscCallBLAS("LAPACKpotri", LAPACKpotri_("L", &B_N, array + cum, &B_N, &B_ierr)); 1825 PetscCheck(!B_ierr, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in POTRI Lapack routine %" PetscBLASInt_FMT, B_ierr); 1826 } else if (use_sytr) { 1827 PetscCallBLAS("LAPACKsytrf", LAPACKsytrf_("L", &B_N, array + cum, &B_N, pivots, Bwork, &B_lwork, &B_ierr)); 1828 PetscCheck(!B_ierr, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in SYTRF Lapack routine %" PetscBLASInt_FMT, B_ierr); 1829 PetscCallBLAS("LAPACKsytri", LAPACKsytri_("L", &B_N, array + cum, &B_N, pivots, Bwork, &B_ierr)); 1830 PetscCheck(!B_ierr, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in SYTRI Lapack routine %" PetscBLASInt_FMT, B_ierr); 1831 } else { 1832 PetscCallBLAS("LAPACKgetrf", LAPACKgetrf_(&B_N, &B_N, array + cum, &B_N, pivots, &B_ierr)); 1833 PetscCheck(!B_ierr, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in GETRF Lapack routine %" PetscBLASInt_FMT, B_ierr); 1834 PetscCallBLAS("LAPACKgetri", LAPACKgetri_(&B_N, array + cum, &B_N, pivots, Bwork, &B_lwork, &B_ierr)); 1835 PetscCheck(!B_ierr, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in GETRI Lapack routine %" PetscBLASInt_FMT, B_ierr); 1836 } 1837 PetscCall(PetscLogFlops(1.0 * subset_size * subset_size * subset_size)); 1838 PetscCall(PetscFPTrapPop()); 1839 cum += subset_size * subset_size; 1840 } 1841 PetscCall(MatSeqAIJRestoreArray(sub_schurs->sum_S_Ej_tilda_all, &array)); 1842 PetscCall(PetscObjectReference((PetscObject)sub_schurs->sum_S_Ej_all)); 1843 PetscCall(MatDestroy(&sub_schurs->sum_S_Ej_inv_all)); 1844 sub_schurs->sum_S_Ej_inv_all = sub_schurs->sum_S_Ej_all; 1845 } 1846 } 1847 PetscCall(VecDestroy(&lstash)); 1848 PetscCall(VecDestroy(&gstash)); 1849 PetscCall(VecScatterDestroy(&sstash)); 1850 1851 if (matl_dbg_viewer) { 1852 if (sub_schurs->S_Ej_all) { 1853 PetscCall(PetscObjectSetName((PetscObject)sub_schurs->S_Ej_all, "SE")); 1854 PetscCall(MatView(sub_schurs->S_Ej_all, matl_dbg_viewer)); 1855 } 1856 if (sub_schurs->sum_S_Ej_all) { 1857 PetscCall(PetscObjectSetName((PetscObject)sub_schurs->sum_S_Ej_all, "SSE")); 1858 PetscCall(MatView(sub_schurs->sum_S_Ej_all, matl_dbg_viewer)); 1859 } 1860 if (sub_schurs->sum_S_Ej_inv_all) { 1861 PetscCall(PetscObjectSetName((PetscObject)sub_schurs->sum_S_Ej_inv_all, "SSEm")); 1862 PetscCall(MatView(sub_schurs->sum_S_Ej_inv_all, matl_dbg_viewer)); 1863 } 1864 if (sub_schurs->sum_S_Ej_tilda_all) { 1865 PetscCall(PetscObjectSetName((PetscObject)sub_schurs->sum_S_Ej_tilda_all, "SSEt")); 1866 PetscCall(MatView(sub_schurs->sum_S_Ej_tilda_all, matl_dbg_viewer)); 1867 } 1868 } 1869 1870 /* when not explicit, we need to set the factor type */ 1871 if (sub_schurs->mat_factor_type == MAT_FACTOR_NONE) sub_schurs->mat_factor_type = sub_schurs->is_hermitian ? MAT_FACTOR_CHOLESKY : MAT_FACTOR_LU; 1872 1873 /* free workspace */ 1874 if (matl_dbg_viewer) PetscCall(PetscViewerFlush(matl_dbg_viewer)); 1875 if (sub_schurs->debug) PetscCallMPI(MPI_Barrier(comm_n)); 1876 PetscCall(PetscViewerDestroy(&matl_dbg_viewer)); 1877 PetscCall(PetscFree2(Bwork, pivots)); 1878 PetscCall(PetscCommDestroy(&comm_n)); 1879 PetscFunctionReturn(PETSC_SUCCESS); 1880 } 1881 1882 PetscErrorCode PCBDDCSubSchursInit(PCBDDCSubSchurs sub_schurs, const char *prefix, IS is_I, IS is_B, PCBDDCGraph graph, ISLocalToGlobalMapping BtoNmap, PetscBool copycc, PetscBool gdsw) 1883 { 1884 IS *faces, *edges, *all_cc, vertices; 1885 PetscInt s, i, n_faces, n_edges, n_all_cc; 1886 PetscBool is_sorted, ispardiso, ismumps; 1887 1888 PetscFunctionBegin; 1889 PetscCall(ISSorted(is_I, &is_sorted)); 1890 PetscCheck(is_sorted, PetscObjectComm((PetscObject)is_I), PETSC_ERR_PLIB, "IS for I dofs should be shorted"); 1891 PetscCall(ISSorted(is_B, &is_sorted)); 1892 PetscCheck(is_sorted, PetscObjectComm((PetscObject)is_B), PETSC_ERR_PLIB, "IS for B dofs should be shorted"); 1893 1894 /* reset any previous data */ 1895 PetscCall(PCBDDCSubSchursReset(sub_schurs)); 1896 1897 sub_schurs->gdsw = gdsw; 1898 1899 /* get index sets for faces and edges (already sorted by global ordering) */ 1900 PetscCall(PCBDDCGraphGetCandidatesIS(graph, &n_faces, &faces, &n_edges, &edges, &vertices)); 1901 n_all_cc = n_faces + n_edges; 1902 PetscCall(PetscBTCreate(n_all_cc, &sub_schurs->is_edge)); 1903 PetscCall(PetscMalloc1(n_all_cc, &all_cc)); 1904 n_all_cc = 0; 1905 for (i = 0; i < n_faces; i++) { 1906 PetscCall(ISGetSize(faces[i], &s)); 1907 if (!s) continue; 1908 if (copycc) { 1909 PetscCall(ISDuplicate(faces[i], &all_cc[n_all_cc])); 1910 } else { 1911 PetscCall(PetscObjectReference((PetscObject)faces[i])); 1912 all_cc[n_all_cc] = faces[i]; 1913 } 1914 n_all_cc++; 1915 } 1916 for (i = 0; i < n_edges; i++) { 1917 PetscCall(ISGetSize(edges[i], &s)); 1918 if (!s) continue; 1919 if (copycc) { 1920 PetscCall(ISDuplicate(edges[i], &all_cc[n_all_cc])); 1921 } else { 1922 PetscCall(PetscObjectReference((PetscObject)edges[i])); 1923 all_cc[n_all_cc] = edges[i]; 1924 } 1925 PetscCall(PetscBTSet(sub_schurs->is_edge, n_all_cc)); 1926 n_all_cc++; 1927 } 1928 PetscCall(PetscObjectReference((PetscObject)vertices)); 1929 sub_schurs->is_vertices = vertices; 1930 PetscCall(PCBDDCGraphRestoreCandidatesIS(graph, &n_faces, &faces, &n_edges, &edges, &vertices)); 1931 sub_schurs->is_dir = NULL; 1932 PetscCall(PCBDDCGraphGetDirichletDofsB(graph, &sub_schurs->is_dir)); 1933 1934 /* Determine if MatFactor can be used */ 1935 PetscCall(PetscStrallocpy(prefix, &sub_schurs->prefix)); 1936 #if defined(PETSC_HAVE_MUMPS) 1937 PetscCall(PetscStrncpy(sub_schurs->mat_solver_type, MATSOLVERMUMPS, sizeof(sub_schurs->mat_solver_type))); 1938 #elif defined(PETSC_HAVE_MKL_PARDISO) 1939 PetscCall(PetscStrncpy(sub_schurs->mat_solver_type, MATSOLVERMKL_PARDISO, sizeof(sub_schurs->mat_solver_type))); 1940 #else 1941 PetscCall(PetscStrncpy(sub_schurs->mat_solver_type, MATSOLVERPETSC, sizeof(sub_schurs->mat_solver_type))); 1942 #endif 1943 sub_schurs->mat_factor_type = MAT_FACTOR_NONE; 1944 #if defined(PETSC_USE_COMPLEX) 1945 sub_schurs->is_hermitian = PETSC_FALSE; /* Hermitian Cholesky is not supported by PETSc and external packages */ 1946 #else 1947 sub_schurs->is_hermitian = PETSC_TRUE; 1948 #endif 1949 sub_schurs->is_posdef = PETSC_TRUE; 1950 sub_schurs->is_symmetric = PETSC_TRUE; 1951 sub_schurs->debug = PETSC_FALSE; 1952 sub_schurs->restrict_comm = PETSC_FALSE; 1953 PetscOptionsBegin(PetscObjectComm((PetscObject)graph->l2gmap), sub_schurs->prefix, "BDDC sub_schurs options", "PC"); 1954 PetscCall(PetscOptionsString("-sub_schurs_mat_solver_type", "Specific direct solver to use", NULL, sub_schurs->mat_solver_type, sub_schurs->mat_solver_type, sizeof(sub_schurs->mat_solver_type), NULL)); 1955 PetscCall(PetscOptionsEnum("-sub_schurs_mat_factor_type", "Factor type to use. Use MAT_FACTOR_NONE for automatic selection", NULL, MatFactorTypes, (PetscEnum)sub_schurs->mat_factor_type, (PetscEnum *)&sub_schurs->mat_factor_type, NULL)); 1956 PetscCall(PetscOptionsBool("-sub_schurs_symmetric", "Symmetric problem", NULL, sub_schurs->is_symmetric, &sub_schurs->is_symmetric, NULL)); 1957 PetscCall(PetscOptionsBool("-sub_schurs_hermitian", "Hermitian problem", NULL, sub_schurs->is_hermitian, &sub_schurs->is_hermitian, NULL)); 1958 PetscCall(PetscOptionsBool("-sub_schurs_posdef", "Positive definite problem", NULL, sub_schurs->is_posdef, &sub_schurs->is_posdef, NULL)); 1959 PetscCall(PetscOptionsBool("-sub_schurs_restrictcomm", "Restrict communicator on active processes only", NULL, sub_schurs->restrict_comm, &sub_schurs->restrict_comm, NULL)); 1960 PetscCall(PetscOptionsBool("-sub_schurs_debug", "Debug output", NULL, sub_schurs->debug, &sub_schurs->debug, NULL)); 1961 PetscOptionsEnd(); 1962 PetscCall(PetscStrcmp(sub_schurs->mat_solver_type, MATSOLVERMUMPS, &ismumps)); 1963 PetscCall(PetscStrcmp(sub_schurs->mat_solver_type, MATSOLVERMKL_PARDISO, &ispardiso)); 1964 sub_schurs->schur_explicit = (PetscBool)(ispardiso || ismumps); 1965 1966 /* for reals, symmetric and Hermitian are synonyms */ 1967 #if !defined(PETSC_USE_COMPLEX) 1968 sub_schurs->is_symmetric = (PetscBool)(sub_schurs->is_symmetric && sub_schurs->is_hermitian); 1969 sub_schurs->is_hermitian = sub_schurs->is_symmetric; 1970 #endif 1971 1972 PetscCall(PetscObjectReference((PetscObject)is_I)); 1973 sub_schurs->is_I = is_I; 1974 PetscCall(PetscObjectReference((PetscObject)is_B)); 1975 sub_schurs->is_B = is_B; 1976 PetscCall(PetscObjectReference((PetscObject)graph->l2gmap)); 1977 sub_schurs->l2gmap = graph->l2gmap; 1978 PetscCall(PetscObjectReference((PetscObject)BtoNmap)); 1979 sub_schurs->BtoNmap = BtoNmap; 1980 sub_schurs->n_subs = n_all_cc; 1981 sub_schurs->is_subs = all_cc; 1982 sub_schurs->S_Ej_all = NULL; 1983 sub_schurs->sum_S_Ej_all = NULL; 1984 sub_schurs->sum_S_Ej_inv_all = NULL; 1985 sub_schurs->sum_S_Ej_tilda_all = NULL; 1986 sub_schurs->is_Ej_all = NULL; 1987 PetscFunctionReturn(PETSC_SUCCESS); 1988 } 1989 1990 PetscErrorCode PCBDDCSubSchursCreate(PCBDDCSubSchurs *sub_schurs) 1991 { 1992 PCBDDCSubSchurs schurs_ctx; 1993 1994 PetscFunctionBegin; 1995 PetscCall(PetscNew(&schurs_ctx)); 1996 schurs_ctx->n_subs = 0; 1997 *sub_schurs = schurs_ctx; 1998 PetscFunctionReturn(PETSC_SUCCESS); 1999 } 2000 2001 PetscErrorCode PCBDDCSubSchursReset(PCBDDCSubSchurs sub_schurs) 2002 { 2003 PetscInt i; 2004 2005 PetscFunctionBegin; 2006 if (!sub_schurs) PetscFunctionReturn(PETSC_SUCCESS); 2007 PetscCall(PetscFree(sub_schurs->prefix)); 2008 PetscCall(MatDestroy(&sub_schurs->A)); 2009 PetscCall(MatDestroy(&sub_schurs->S)); 2010 PetscCall(ISDestroy(&sub_schurs->is_I)); 2011 PetscCall(ISDestroy(&sub_schurs->is_B)); 2012 PetscCall(ISLocalToGlobalMappingDestroy(&sub_schurs->l2gmap)); 2013 PetscCall(ISLocalToGlobalMappingDestroy(&sub_schurs->BtoNmap)); 2014 PetscCall(MatDestroy(&sub_schurs->S_Ej_all)); 2015 PetscCall(MatDestroy(&sub_schurs->sum_S_Ej_all)); 2016 PetscCall(MatDestroy(&sub_schurs->sum_S_Ej_inv_all)); 2017 PetscCall(MatDestroy(&sub_schurs->sum_S_Ej_tilda_all)); 2018 PetscCall(ISDestroy(&sub_schurs->is_Ej_all)); 2019 PetscCall(ISDestroy(&sub_schurs->is_vertices)); 2020 PetscCall(ISDestroy(&sub_schurs->is_dir)); 2021 PetscCall(PetscBTDestroy(&sub_schurs->is_edge)); 2022 for (i = 0; i < sub_schurs->n_subs; i++) PetscCall(ISDestroy(&sub_schurs->is_subs[i])); 2023 if (sub_schurs->n_subs) PetscCall(PetscFree(sub_schurs->is_subs)); 2024 if (sub_schurs->reuse_solver) PetscCall(PCBDDCReuseSolversReset(sub_schurs->reuse_solver)); 2025 PetscCall(PetscFree(sub_schurs->reuse_solver)); 2026 if (sub_schurs->change) { 2027 for (i = 0; i < sub_schurs->n_subs; i++) { 2028 PetscCall(KSPDestroy(&sub_schurs->change[i])); 2029 PetscCall(ISDestroy(&sub_schurs->change_primal_sub[i])); 2030 } 2031 } 2032 PetscCall(PetscFree(sub_schurs->change)); 2033 PetscCall(PetscFree(sub_schurs->change_primal_sub)); 2034 sub_schurs->n_subs = 0; 2035 PetscFunctionReturn(PETSC_SUCCESS); 2036 } 2037 2038 PetscErrorCode PCBDDCSubSchursDestroy(PCBDDCSubSchurs *sub_schurs) 2039 { 2040 PetscFunctionBegin; 2041 PetscCall(PCBDDCSubSchursReset(*sub_schurs)); 2042 PetscCall(PetscFree(*sub_schurs)); 2043 PetscFunctionReturn(PETSC_SUCCESS); 2044 } 2045 2046 static inline PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt *queue_tip, PetscInt n_prev, PetscBT touched, PetscInt *xadj, PetscInt *adjncy, PetscInt *n_added) 2047 { 2048 PetscInt i, j, n; 2049 2050 PetscFunctionBegin; 2051 n = 0; 2052 for (i = -n_prev; i < 0; i++) { 2053 PetscInt start_dof = queue_tip[i]; 2054 for (j = xadj[start_dof]; j < xadj[start_dof + 1]; j++) { 2055 PetscInt dof = adjncy[j]; 2056 if (!PetscBTLookup(touched, dof)) { 2057 PetscCall(PetscBTSet(touched, dof)); 2058 queue_tip[n] = dof; 2059 n++; 2060 } 2061 } 2062 } 2063 *n_added = n; 2064 PetscFunctionReturn(PETSC_SUCCESS); 2065 } 2066